Package madgraph :: Package madevent :: Module gen_ximprove
[hide private]
[frames] | no frames]

Source Code for Module madgraph.madevent.gen_ximprove

   1  ################################################################################ 
   2  # 
   3  # Copyright (c) 2014 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  """ A python file to replace the fortran script gen_ximprove. 
  16      This script analyses the result of the survey/ previous refine and  
  17      creates the jobs for the following script. 
  18  """ 
  19  from __future__ import division 
  20   
  21  import collections 
  22  import os 
  23  import glob 
  24  import logging 
  25  import math 
  26  import re 
  27  import subprocess 
  28  import shutil 
  29  import stat 
  30  import sys 
  31   
  32  try: 
  33      import madgraph 
  34  except ImportError: 
  35      MADEVENT = True 
  36      import internal.sum_html as sum_html 
  37      import internal.banner as bannermod 
  38      import internal.misc as misc 
  39      import internal.files as files 
  40      import internal.cluster as cluster 
  41      import internal.combine_grid as combine_grid 
  42      import internal.combine_runs as combine_runs 
  43      import internal.lhe_parser as lhe_parser 
  44  else: 
  45      MADEVENT= False 
  46      import madgraph.madevent.sum_html as sum_html 
  47      import madgraph.various.banner as bannermod 
  48      import madgraph.various.misc as misc 
  49      import madgraph.iolibs.files as files 
  50      import madgraph.various.cluster as cluster 
  51      import madgraph.madevent.combine_grid as combine_grid 
  52      import madgraph.madevent.combine_runs as combine_runs 
  53      import madgraph.various.lhe_parser as lhe_parser 
  54   
  55  logger = logging.getLogger('madgraph.madevent.gen_ximprove') 
  56  pjoin = os.path.join 
57 58 -class gensym(object):
59 """a class to call the fortran gensym executable and handle it's output 60 in order to create the various job that are needed for the survey""" 61 62 #convenient shortcut for the formatting of variable 63 @ staticmethod
64 - def format_variable(*args):
65 return bannermod.ConfigFile.format_variable(*args)
66 67 combining_job = 2 # number of channel by ajob 68 splitted_grid = False 69 min_iterations = 3 70 mode= "survey" 71 72
73 - def __init__(self, cmd, opt=None):
74 75 try: 76 super(gensym, self).__init__(cmd, opt) 77 except TypeError: 78 pass 79 80 # Run statistics, a dictionary of RunStatistics(), with 81 self.run_statistics = {} 82 83 self.cmd = cmd 84 self.run_card = cmd.run_card 85 self.me_dir = cmd.me_dir 86 87 88 # dictionary to keep track of the precision when combining iteration 89 self.cross = collections.defaultdict(int) 90 self.abscross = collections.defaultdict(int) 91 self.sigma = collections.defaultdict(int) 92 self.chi2 = collections.defaultdict(int) 93 94 self.splitted_grid = False 95 if self.cmd.proc_characteristics['loop_induced']: 96 nexternal = self.cmd.proc_characteristics['nexternal'] 97 self.splitted_grid = max(2, (nexternal-2)**2) 98 if hasattr(self.cmd, "opts") and self.cmd.opts['accuracy'] == 0.1: 99 self.cmd.opts['accuracy'] = 0.02 100 101 if isinstance(cmd.cluster, cluster.MultiCore) and self.splitted_grid > 1: 102 self.splitted_grid = int(cmd.cluster.nb_core**0.5) 103 if self.splitted_grid == 1 and cmd.cluster.nb_core >1: 104 self.splitted_grid = 2 105 106 #if the user defines it in the run_card: 107 if self.run_card['survey_splitting'] != -1: 108 self.splitted_grid = self.run_card['survey_splitting'] 109 if self.run_card['survey_nchannel_per_job'] != -1: 110 self.combining_job = self.run_card['survey_nchannel_per_job'] 111 112 self.splitted_Pdir = {} 113 self.splitted_for_dir = lambda x,y: self.splitted_grid 114 self.combining_job_for_Pdir = lambda x: self.combining_job 115 self.lastoffset = {}
116
117 - def launch(self, to_submit=True, clean=True):
118 """ """ 119 120 self.subproc = [l.strip() for l in open(pjoin(self.me_dir,'SubProcesses', 121 'subproc.mg'))] 122 subproc = self.subproc 123 124 P_zero_result = [] # check the number of times where they are no phase-space 125 126 nb_tot_proc = len(subproc) 127 job_list = {} 128 for nb_proc,subdir in enumerate(subproc): 129 self.cmd.update_status('Compiling for process %s/%s. <br> (previous processes already running)' % \ 130 (nb_proc+1,nb_tot_proc), level=None) 131 132 subdir = subdir.strip() 133 Pdir = pjoin(self.me_dir, 'SubProcesses',subdir) 134 logger.info(' %s ' % subdir) 135 136 # clean previous run 137 if clean: 138 for match in misc.glob('*ajob*', Pdir): 139 if os.path.basename(match)[:4] in ['ajob', 'wait', 'run.', 'done']: 140 os.remove(match) 141 for match in misc.glob('G*', Pdir): 142 if os.path.exists(pjoin(match,'results.dat')): 143 os.remove(pjoin(match, 'results.dat')) 144 if os.path.exists(pjoin(match, 'ftn25')): 145 os.remove(pjoin(match, 'ftn25')) 146 147 #compile gensym 148 self.cmd.compile(['gensym'], cwd=Pdir) 149 if not os.path.exists(pjoin(Pdir, 'gensym')): 150 raise Exception, 'Error make gensym not successful' 151 152 # Launch gensym 153 p = misc.Popen(['./gensym'], stdout=subprocess.PIPE, 154 stderr=subprocess.STDOUT, cwd=Pdir) 155 #sym_input = "%(points)d %(iterations)d %(accuracy)f \n" % self.opts 156 (stdout, _) = p.communicate('') 157 158 if os.path.exists(pjoin(self.me_dir,'error')): 159 files.mv(pjoin(self.me_dir,'error'), pjoin(Pdir,'ajob.no_ps.log')) 160 P_zero_result.append(subdir) 161 continue 162 163 jobs = stdout.split() 164 job_list[Pdir] = jobs 165 try: 166 # check that all input are valid 167 [float(s) for s in jobs] 168 except Exception: 169 logger.debug("unformated string found in gensym. Please check:\n %s" % stdout) 170 done=False 171 job_list[Pdir] = [] 172 lines = stdout.split('\n') 173 for l in lines: 174 try: 175 [float(s) for s in l.split()] 176 except: 177 continue 178 else: 179 if done: 180 raise Exception, 'Parsing error in gensym: %s' % stdout 181 job_list[Pdir] = l.split() 182 done = True 183 if not done: 184 raise Exception, 'Parsing error in gensym: %s' % stdout 185 186 self.cmd.compile(['madevent'], cwd=Pdir) 187 if to_submit: 188 self.submit_to_cluster(job_list) 189 job_list = {} 190 191 return job_list, P_zero_result
192
193 - def resubmit(self, min_precision=1.0, resubmit_zero=False):
194 """collect the result of the current run and relaunch each channel 195 not completed or optionally a completed one with a precision worse than 196 a threshold (and/or the zero result channel)""" 197 198 job_list, P_zero_result = self.launch(to_submit=False, clean=False) 199 200 for P , jobs in dict(job_list).items(): 201 misc.sprint(jobs) 202 to_resub = [] 203 for job in jobs: 204 if os.path.exists(pjoin(P, 'G%s' % job)) and os.path.exists(pjoin(P, 'G%s' % job, 'results.dat')): 205 one_result = sum_html.OneResult(job) 206 try: 207 one_result.read_results(pjoin(P, 'G%s' % job, 'results.dat')) 208 except: 209 to_resub.append(job) 210 if one_result.xsec == 0: 211 if resubmit_zero: 212 to_resub.append(job) 213 elif max(one_result.xerru, one_result.xerrc)/one_result.xsec > min_precision: 214 to_resub.append(job) 215 else: 216 to_resub.append(job) 217 if to_resub: 218 for G in to_resub: 219 try: 220 shutil.rmtree(pjoin(P, 'G%s' % G)) 221 except Exception, error: 222 misc.sprint(error) 223 pass 224 misc.sprint(to_resub) 225 self.submit_to_cluster({P: to_resub})
226 227 228 229 230 231 232 233 234 235 236
237 - def submit_to_cluster(self, job_list):
238 """ """ 239 240 if self.run_card['job_strategy'] > 0: 241 if len(job_list) >1: 242 for path, dirs in job_list.items(): 243 self.submit_to_cluster({path:dirs}) 244 return 245 path, value = job_list.items()[0] 246 nexternal = self.cmd.proc_characteristics['nexternal'] 247 current = open(pjoin(path, "nexternal.inc")).read() 248 ext = re.search(r"PARAMETER \(NEXTERNAL=(\d+)\)", current).group(1) 249 250 if self.run_card['job_strategy'] == 2: 251 self.splitted_grid = 2 252 if nexternal == int(ext): 253 to_split = 2 254 else: 255 to_split = 0 256 if hasattr(self, 'splitted_Pdir'): 257 self.splitted_Pdir[path] = to_split 258 else: 259 self.splitted_Pdir = {path: to_split} 260 self.splitted_for_dir = lambda x,y : self.splitted_Pdir[x] 261 elif self.run_card['job_strategy'] == 1: 262 if nexternal == int(ext): 263 combine = 1 264 else: 265 combine = self.combining_job 266 if hasattr(self, 'splitted_Pdir'): 267 self.splitted_Pdir[path] = combine 268 else: 269 self.splitted_Pdir = {path: combine} 270 self.combining_job_for_Pdir = lambda x : self.splitted_Pdir[x] 271 272 if not self.splitted_grid: 273 return self.submit_to_cluster_no_splitting(job_list) 274 elif self.cmd.cluster_mode == 0: 275 return self.submit_to_cluster_no_splitting(job_list) 276 elif self.cmd.cluster_mode == 2 and self.cmd.options['nb_core'] == 1: 277 return self.submit_to_cluster_no_splitting(job_list) 278 else: 279 return self.submit_to_cluster_splitted(job_list)
280 281
282 - def submit_to_cluster_no_splitting(self, job_list):
283 """submit the survey without the parralelization. 284 This is the old mode which is still usefull in single core""" 285 286 # write the template file for the parameter file 287 self.write_parameter(parralelization=False, Pdirs=job_list.keys()) 288 289 290 # launch the job with the appropriate grouping 291 for Pdir, jobs in job_list.items(): 292 jobs = list(jobs) 293 i=0 294 while jobs: 295 i+=1 296 to_submit = ['0'] # the first entry is actually the offset 297 for _ in range(self.combining_job_for_Pdir(Pdir)): 298 if jobs: 299 to_submit.append(jobs.pop(0)) 300 301 self.cmd.launch_job(pjoin(self.me_dir, 'SubProcesses', 'survey.sh'), 302 argument=to_submit, 303 cwd=pjoin(self.me_dir,'SubProcesses' , Pdir))
304 305
306 - def create_resubmit_one_iter(self, Pdir, G, submit_ps, nb_job, step=0):
307 """prepare the input_file for submitting the channel""" 308 309 310 if 'SubProcesses' not in Pdir: 311 Pdir = pjoin(self.me_dir, 'SubProcesses', Pdir) 312 313 #keep track of how many job are sended 314 self.splitted_Pdir[(Pdir, G)] = int(nb_job) 315 316 317 # 1. write the new input_app.txt 318 run_card = self.cmd.run_card 319 options = {'event' : submit_ps, 320 'maxiter': 1, 321 'miniter': 1, 322 'accuracy': self.cmd.opts['accuracy'], 323 'helicity': run_card['nhel_survey'] if 'nhel_survey' in run_card \ 324 else run_card['nhel'], 325 'gridmode': -2, 326 'channel' : G 327 } 328 329 Gdir = pjoin(Pdir, 'G%s' % G) 330 self.write_parameter_file(pjoin(Gdir, 'input_app.txt'), options) 331 332 # 2. check that ftn25 exists. 333 assert os.path.exists(pjoin(Gdir, "ftn25")) 334 335 336 # 3. Submit the new jobs 337 #call back function 338 packet = cluster.Packet((Pdir, G, step+1), 339 self.combine_iteration, 340 (Pdir, G, step+1)) 341 342 if step ==0: 343 self.lastoffset[(Pdir, G)] = 0 344 345 # resubmit the new jobs 346 for i in xrange(int(nb_job)): 347 name = "G%s_%s" % (G,i+1) 348 self.lastoffset[(Pdir, G)] += 1 349 offset = self.lastoffset[(Pdir, G)] 350 self.cmd.launch_job(pjoin(self.me_dir, 'SubProcesses', 'refine_splitted.sh'), 351 argument=[name, 'G%s'%G, offset], 352 cwd= Pdir, 353 packet_member=packet)
354 355
356 - def submit_to_cluster_splitted(self, job_list):
357 """ submit the version of the survey with splitted grid creation 358 """ 359 360 #if self.splitted_grid <= 1: 361 # return self.submit_to_cluster_no_splitting(job_list) 362 363 for Pdir, jobs in job_list.items(): 364 if not jobs: 365 continue 366 if self.splitted_for_dir(Pdir, jobs[0]) <= 1: 367 return self.submit_to_cluster_no_splitting({Pdir:jobs}) 368 369 self.write_parameter(parralelization=True, Pdirs=[Pdir]) 370 # launch the job with the appropriate grouping 371 372 for job in jobs: 373 packet = cluster.Packet((Pdir, job, 1), self.combine_iteration, (Pdir, job, 1)) 374 for i in range(self.splitted_for_dir(Pdir, job)): 375 self.cmd.launch_job(pjoin(self.me_dir, 'SubProcesses', 'survey.sh'), 376 argument=[i+1, job], 377 cwd=pjoin(self.me_dir,'SubProcesses' , Pdir), 378 packet_member=packet)
379
380 - def combine_iteration(self, Pdir, G, step):
381 382 grid_calculator, cross, error = self.combine_grid(Pdir, G, step) 383 384 # Compute the number of events used for this run. 385 nb_events = grid_calculator.target_evt 386 387 Gdirs = [] #build the the list of directory 388 for i in range(self.splitted_for_dir(Pdir, G)): 389 path = pjoin(Pdir, "G%s_%s" % (G, i+1)) 390 Gdirs.append(path) 391 392 # 4. make the submission of the next iteration 393 # Three cases - less than 3 iteration -> continue 394 # - more than 3 and less than 5 -> check error 395 # - more than 5 -> prepare info for refine 396 need_submit = False 397 if step < self.min_iterations and cross != 0: 398 if step == 1: 399 need_submit = True 400 else: 401 across = self.abscross[(Pdir,G)]/(self.sigma[(Pdir,G)]+1e-99) 402 tot_across = self.get_current_axsec() 403 if across / tot_across < 1e-6: 404 need_submit = False 405 elif error < self.cmd.opts['accuracy'] / 100: 406 need_submit = False 407 else: 408 need_submit = True 409 410 elif step >= self.cmd.opts['iterations']: 411 need_submit = False 412 elif self.cmd.opts['accuracy'] < 0: 413 #check for luminosity 414 raise Exception, "Not Implemented" 415 elif self.abscross[(Pdir,G)] == 0: 416 need_submit = False 417 else: 418 across = self.abscross[(Pdir,G)]/(self.sigma[(Pdir,G)]+1e-99) 419 tot_across = self.get_current_axsec() 420 if across == 0: 421 need_submit = False 422 elif across / tot_across < 1e-5: 423 need_submit = False 424 elif error > self.cmd.opts['accuracy']: 425 need_submit = True 426 else: 427 need_submit = False 428 429 430 if cross: 431 grid_calculator.write_grid_for_submission(Pdir,G, 432 self.splitted_for_dir(Pdir, G), 433 nb_events,mode=self.mode, 434 conservative_factor=5.0) 435 436 xsec_format = '.%ig'%(max(3,int(math.log10(1.0/float(error)))+2) 437 if float(cross)!=0.0 and float(error)!=0.0 else 8) 438 if need_submit: 439 message = "%%s/G%%s is at %%%s +- %%.3g pb. Now submitting iteration #%s."%(xsec_format, step+1) 440 logger.info(message%\ 441 (os.path.basename(Pdir), G, float(cross), 442 float(error)*float(cross))) 443 self.resubmit_survey(Pdir,G, Gdirs, step) 444 elif cross: 445 logger.info("Survey finished for %s/G%s at %s"%( 446 os.path.basename(Pdir),G,('%%%s +- %%.3g pb'%xsec_format))% 447 (float(cross), float(error)*float(cross))) 448 # prepare information for refine 449 newGpath = pjoin(self.me_dir,'SubProcesses' , Pdir, 'G%s' % G) 450 if not os.path.exists(newGpath): 451 os.mkdir(newGpath) 452 453 # copy the new grid: 454 files.cp(pjoin(Gdirs[0], 'ftn25'), 455 pjoin(self.me_dir,'SubProcesses' , Pdir, 'G%s' % G, 'ftn26')) 456 457 # copy the events 458 fsock = open(pjoin(newGpath, 'events.lhe'), 'w') 459 for Gdir in Gdirs: 460 fsock.write(open(pjoin(Gdir, 'events.lhe')).read()) 461 462 # copy one log 463 files.cp(pjoin(Gdirs[0], 'log.txt'), 464 pjoin(self.me_dir,'SubProcesses' , Pdir, 'G%s' % G)) 465 466 467 # create the appropriate results.dat 468 self.write_results(grid_calculator, cross, error, Pdir, G, step) 469 else: 470 logger.info("Survey finished for %s/G%s [0 cross]", os.path.basename(Pdir),G) 471 472 Gdir = pjoin(self.me_dir,'SubProcesses' , Pdir, 'G%s' % G) 473 if not os.path.exists(Gdir): 474 os.mkdir(Gdir) 475 # copy one log 476 files.cp(pjoin(Gdirs[0], 'log.txt'), Gdir) 477 # create the appropriate results.dat 478 self.write_results(grid_calculator, cross, error, Pdir, G, step) 479 480 return 0
481
482 - def combine_grid(self, Pdir, G, step, exclude_sub_jobs=[]):
483 """ exclude_sub_jobs is to remove some of the subjobs if a numerical 484 issue is detected in one of them. Warning is issue when this occurs. 485 """ 486 487 # 1. create an object to combine the grid information and fill it 488 grid_calculator = combine_grid.grid_information(self.run_card['nhel']) 489 490 for i in range(self.splitted_for_dir(Pdir, G)): 491 if i in exclude_sub_jobs: 492 continue 493 path = pjoin(Pdir, "G%s_%s" % (G, i+1)) 494 fsock = misc.mult_try_open(pjoin(path, 'results.dat')) 495 one_result = grid_calculator.add_results_information(fsock) 496 fsock.close() 497 if one_result.axsec == 0: 498 grid_calculator.onefail = True 499 continue # grid_information might not exists 500 fsock = misc.mult_try_open(pjoin(path, 'grid_information')) 501 grid_calculator.add_one_grid_information(fsock) 502 fsock.close() 503 os.remove(pjoin(path, 'results.dat')) 504 #os.remove(pjoin(path, 'grid_information')) 505 506 507 508 #2. combine the information about the total crossection / error 509 # start by keep the interation in memory 510 cross, across, sigma = grid_calculator.get_cross_section() 511 512 #3. Try to avoid one single PS point which ruins the integration 513 # Should be related to loop evaluation instability. 514 maxwgt = grid_calculator.get_max_wgt(0.01) 515 if maxwgt: 516 nunwgt = grid_calculator.get_nunwgt(maxwgt) 517 # Make sure not to apply the security below during the first step of the 518 # survey. Also, disregard channels with a contribution relative to the 519 # total cross-section smaller than 1e-8 since in this case it is unlikely 520 # that this channel will need more than 1 event anyway. 521 apply_instability_security = False 522 rel_contrib = 0.0 523 if (self.__class__ != gensym or step > 1): 524 Pdir_across = 0.0 525 Gdir_across = 0.0 526 for (mPdir,mG) in self.abscross.keys(): 527 if mPdir == Pdir: 528 Pdir_across += (self.abscross[(mPdir,mG)]/ 529 (self.sigma[(mPdir,mG)]+1e-99)) 530 if mG == G: 531 Gdir_across += (self.abscross[(mPdir,mG)]/ 532 (self.sigma[(mPdir,mG)]+1e-99)) 533 rel_contrib = abs(Gdir_across/(Pdir_across+1e-99)) 534 if rel_contrib > (1.0e-8) and \ 535 nunwgt < 2 and len(grid_calculator.results) > 1: 536 apply_instability_security = True 537 538 if apply_instability_security: 539 # check the ratio between the different submit 540 th_maxwgt = [(r.th_maxwgt,i) for i,r in enumerate(grid_calculator.results)] 541 th_maxwgt.sort() 542 ratio = th_maxwgt[-1][0]/th_maxwgt[-2][0] 543 if ratio > 1e4: 544 logger.warning( 545 """"One Event with large weight have been found (ratio = %.3g) in channel G%s (with rel.contrib=%.3g). 546 This is likely due to numerical instabilities. The associated job is discarded to recover. 547 For offline investigation, the problematic discarded events are stored in: 548 %s"""%(ratio,G,rel_contrib,pjoin(Pdir,'DiscardedUnstableEvents'))) 549 exclude_sub_jobs = list(exclude_sub_jobs) 550 exclude_sub_jobs.append(th_maxwgt[-1][1]) 551 grid_calculator.results.run_statistics['skipped_subchannel'] += 1 552 553 # Add some monitoring of the problematic events 554 gPath = pjoin(Pdir, "G%s_%s" % (G, th_maxwgt[-1][1]+1)) 555 if os.path.isfile(pjoin(gPath,'events.lhe')): 556 lhe_file = lhe_parser.EventFile(pjoin(gPath,'events.lhe')) 557 discardedPath = pjoin(Pdir,'DiscardedUnstableEvents') 558 if not os.path.exists(discardedPath): 559 os.mkdir(discardedPath) 560 if os.path.isdir(discardedPath): 561 # Keep only the event with a maximum weight, as it surely 562 # is the problematic one. 563 evtRecord = open(pjoin(discardedPath,'discarded_G%s.dat'%G),'a') 564 lhe_file.seek(0) #rewind the file 565 try: 566 evtRecord.write('\n'+str(max(lhe_file,key=lambda evt:abs(evt.wgt)))) 567 except Exception: 568 #something wrong write the full file. 569 lhe_file.close() 570 evtRecord.write(pjoin(gPath,'events.lhe').read()) 571 evtRecord.close() 572 573 return self.combine_grid(Pdir, G, step, exclude_sub_jobs) 574 575 576 if across !=0: 577 if sigma != 0: 578 self.cross[(Pdir,G)] += cross**3/sigma**2 579 self.abscross[(Pdir,G)] += across * cross**2/sigma**2 580 self.sigma[(Pdir,G)] += cross**2/ sigma**2 581 self.chi2[(Pdir,G)] += cross**4/sigma**2 582 # and use those iteration to get the current estimator 583 cross = self.cross[(Pdir,G)]/self.sigma[(Pdir,G)] 584 if step > 1: 585 error = math.sqrt(abs((self.chi2[(Pdir,G)]/cross**2 - \ 586 self.sigma[(Pdir,G)])/(step-1))/self.sigma[(Pdir,G)]) 587 else: 588 error = sigma/cross 589 else: 590 self.cross[(Pdir,G)] = cross 591 self.abscross[(Pdir,G)] = across 592 self.sigma[(Pdir,G)] = 0 593 self.chi2[(Pdir,G)] = 0 594 cross = self.cross[(Pdir,G)] 595 error = 0 596 597 else: 598 error = 0 599 600 grid_calculator.results.compute_values(update_statistics=True) 601 if (str(os.path.basename(Pdir)), G) in self.run_statistics: 602 self.run_statistics[(str(os.path.basename(Pdir)), G)]\ 603 .aggregate_statistics(grid_calculator.results.run_statistics) 604 else: 605 self.run_statistics[(str(os.path.basename(Pdir)), G)] = \ 606 grid_calculator.results.run_statistics 607 608 self.warnings_from_statistics(G, grid_calculator.results.run_statistics) 609 stats_msg = grid_calculator.results.run_statistics.nice_output( 610 '/'.join([os.path.basename(Pdir),'G%s'%G])) 611 612 if stats_msg: 613 logger.log(5, stats_msg) 614 615 # Clean up grid_information to avoid border effects in case of a crash 616 for i in range(self.splitted_for_dir(Pdir, G)): 617 path = pjoin(Pdir, "G%s_%s" % (G, i+1)) 618 try: 619 os.remove(pjoin(path, 'grid_information')) 620 except OSError, oneerror: 621 if oneerror.errno != 2: 622 raise 623 return grid_calculator, cross, error
624
625 - def warnings_from_statistics(self,G,stats):
626 """Possible warn user for worrying MadLoop stats for this channel.""" 627 628 if stats['n_madloop_calls']==0: 629 return 630 631 EPS_fraction = float(stats['exceptional_points'])/stats['n_madloop_calls'] 632 633 msg = "Channel %s has encountered a fraction of %.3g\n"+ \ 634 "of numerically unstable loop matrix element computations\n"+\ 635 "(which could not be rescued using quadruple precision).\n"+\ 636 "The results might not be trusted." 637 638 if 0.01 > EPS_fraction > 0.001: 639 logger.warning(msg%(G,EPS_fraction)) 640 elif EPS_fraction > 0.01: 641 logger.critical((msg%(G,EPS_fraction)).replace('might', 'can')) 642 raise Exception, (msg%(G,EPS_fraction)).replace('might', 'can')
643
644 - def get_current_axsec(self):
645 646 across = 0 647 for (Pdir,G) in self.abscross: 648 across += self.abscross[(Pdir,G)]/(self.sigma[(Pdir,G)]+1e-99) 649 return across
650
651 - def write_results(self, grid_calculator, cross, error, Pdir, G, step):
652 653 #compute the value 654 if cross == 0: 655 abscross,nw, luminosity = 0, 0, 0 656 wgt, maxit,nunwgt, wgt, nevents = 0,0,0,0,0 657 maxwgt = 0 658 error = 0 659 else: 660 grid_calculator.results.compute_values() 661 abscross = self.abscross[(Pdir,G)]/self.sigma[(Pdir,G)] 662 nw = grid_calculator.results.nw 663 wgt = grid_calculator.results.wgt 664 maxit = step 665 wgt = 0 666 nevents = grid_calculator.results.nevents 667 maxwgt = grid_calculator.get_max_wgt() 668 nunwgt = grid_calculator.get_nunwgt() 669 luminosity = nunwgt/cross 670 671 #format the results.dat 672 def fstr(nb): 673 data = '%E' % nb 674 nb, power = data.split('E') 675 nb = float(nb) /10 676 power = int(power) + 1 677 return '%.5fE%+03i' %(nb,power)
678 line = '%s %s %s %i %i %i %i %s %s %s %s 0.0 0\n' % \ 679 (fstr(cross), fstr(error*cross), fstr(error*cross), 680 nevents, nw, maxit,nunwgt, 681 fstr(luminosity), fstr(wgt), fstr(abscross), fstr(maxwgt)) 682 683 fsock = open(pjoin(self.me_dir,'SubProcesses' , Pdir, 'G%s' % G, 684 'results.dat'),'w') 685 fsock.writelines(line) 686 fsock.close()
687
688 - def resubmit_survey(self, Pdir, G, Gdirs, step):
689 """submit the next iteration of the survey""" 690 691 # 1. write the new input_app.txt to double the number of points 692 run_card = self.cmd.run_card 693 options = {'event' : 2**(step) * self.cmd.opts['points'] / self.splitted_grid, 694 'maxiter': 1, 695 'miniter': 1, 696 'accuracy': self.cmd.opts['accuracy'], 697 'helicity': run_card['nhel_survey'] if 'nhel_survey' in run_card \ 698 else run_card['nhel'], 699 'gridmode': -2, 700 'channel' : '' 701 } 702 703 if int(options['helicity']) == 1: 704 options['event'] = options['event'] * 2**(self.cmd.proc_characteristics['nexternal']//3) 705 706 for Gdir in Gdirs: 707 self.write_parameter_file(pjoin(Gdir, 'input_app.txt'), options) 708 709 710 #2. resubmit the new jobs 711 packet = cluster.Packet((Pdir, G, step+1), self.combine_iteration, \ 712 (Pdir, G, step+1)) 713 nb_step = len(Gdirs) * (step+1) 714 for i,subdir in enumerate(Gdirs): 715 subdir = subdir.rsplit('_',1)[1] 716 subdir = int(subdir) 717 offset = nb_step+i+1 718 offset=str(offset) 719 tag = "%s.%s" % (subdir, offset) 720 721 self.cmd.launch_job(pjoin(self.me_dir, 'SubProcesses', 'survey.sh'), 722 argument=[tag, G], 723 cwd=pjoin(self.me_dir,'SubProcesses' , Pdir), 724 packet_member=packet)
725 726 727 728
729 - def write_parameter_file(self, path, options):
730 """ """ 731 732 template =""" %(event)s %(maxiter)s %(miniter)s !Number of events and max and min iterations 733 %(accuracy)s !Accuracy 734 %(gridmode)s !Grid Adjustment 0=none, 2=adjust 735 1 !Suppress Amplitude 1=yes 736 %(helicity)s !Helicity Sum/event 0=exact 737 %(channel)s """ 738 options['event'] = int(options['event']) 739 open(path, 'w').write(template % options)
740 741 742
743 - def write_parameter(self, parralelization, Pdirs=None):
744 """Write the parameter of the survey run""" 745 746 run_card = self.cmd.run_card 747 748 options = {'event' : self.cmd.opts['points'], 749 'maxiter': self.cmd.opts['iterations'], 750 'miniter': self.min_iterations, 751 'accuracy': self.cmd.opts['accuracy'], 752 'helicity': run_card['nhel_survey'] if 'nhel_survey' in run_card \ 753 else run_card['nhel'], 754 'gridmode': 2, 755 'channel': '' 756 } 757 758 if int(options['helicity'])== 1: 759 options['event'] = options['event'] * 2**(self.cmd.proc_characteristics['nexternal']//3) 760 761 if parralelization: 762 options['gridmode'] = -2 763 options['maxiter'] = 1 #this is automatic in dsample anyway 764 options['miniter'] = 1 #this is automatic in dsample anyway 765 options['event'] /= self.splitted_grid 766 767 if not Pdirs: 768 Pdirs = self.subproc 769 770 for Pdir in Pdirs: 771 path =pjoin(Pdir, 'input_app.txt') 772 self.write_parameter_file(path, options)
773
774 775 776 -class gen_ximprove(object):
777 778 779 # some hardcoded value which impact the generation 780 gen_events_security = 1.2 # multiply the number of requested event by this number for security 781 combining_job = 0 # allow to run multiple channel in sequence 782 max_request_event = 1000 # split jobs if a channel if it needs more than that 783 max_event_in_iter = 5000 784 min_event_in_iter = 1000 785 max_splitting = 130 # maximum duplication of a given channel 786 min_iter = 3 787 max_iter = 9 788 keep_grid_for_refine = False # only apply if needed to split the job 789 790 #convenient shortcut for the formatting of variable 791 @ staticmethod
792 - def format_variable(*args):
793 return bannermod.ConfigFile.format_variable(*args)
794 795
796 - def __new__(cls, cmd, opt):
797 """Choose in which type of refine we want to be""" 798 799 if cmd.proc_characteristics['loop_induced']: 800 return super(gen_ximprove, cls).__new__(gen_ximprove_share, cmd, opt) 801 elif gen_ximprove.format_variable(cmd.run_card['gridpack'], bool): 802 return super(gen_ximprove, cls).__new__(gen_ximprove_gridpack, cmd, opt) 803 elif cmd.run_card["job_strategy"] == 2: 804 return super(gen_ximprove, cls).__new__(gen_ximprove_share, cmd, opt) 805 else: 806 return super(gen_ximprove, cls).__new__(gen_ximprove_v4, cmd, opt)
807 808
809 - def __init__(self, cmd, opt=None):
810 811 try: 812 super(gen_ximprove, self).__init__(cmd, opt) 813 except TypeError: 814 pass 815 816 self.run_statistics = {} 817 self.cmd = cmd 818 self.run_card = cmd.run_card 819 run_card = self.run_card 820 self.me_dir = cmd.me_dir 821 822 #extract from the run_card the information that we need. 823 self.gridpack = run_card['gridpack'] 824 self.nhel = run_card['nhel'] 825 if "nhel_refine" in run_card: 826 self.nhel = run_card["nhel_refine"] 827 828 if self.run_card['refine_evt_by_job'] != -1: 829 self.max_request_event = run_card['refine_evt_by_job'] 830 831 832 # Default option for the run 833 self.gen_events = True 834 self.parralel = False 835 # parameter which was input for the normal gen_ximprove run 836 self.err_goal = 0.01 837 self.max_np = 9 838 self.split_channels = False 839 # parameter for the gridpack run 840 self.nreq = 2000 841 self.iseed = 4321 842 843 # placeholder for information 844 self.results = 0 #updated in launch/update_html 845 846 if isinstance(opt, dict): 847 self.configure(opt) 848 elif isinstance(opt, bannermod.GridpackCard): 849 self.configure_gridpack(opt)
850
851 - def __call__(self):
852 return self.launch()
853
854 - def launch(self):
855 """running """ 856 857 #start the run 858 self.handle_seed() 859 self.results = sum_html.collect_result(self.cmd, 860 main_dir=pjoin(self.cmd.me_dir,'SubProcesses')) #main_dir is for gridpack readonly mode 861 if self.gen_events: 862 # We run to provide a given number of events 863 self.get_job_for_event() 864 else: 865 # We run to achieve a given precision 866 self.get_job_for_precision()
867 868
869 - def configure(self, opt):
870 """Defines some parameter of the run""" 871 872 for key, value in opt.items(): 873 if key in self.__dict__: 874 targettype = type(getattr(self, key)) 875 setattr(self, key, self.format_variable(value, targettype, key)) 876 else: 877 raise Exception, '%s not define' % key 878 879 880 # special treatment always do outside the loop to avoid side effect 881 if 'err_goal' in opt: 882 if self.err_goal < 1: 883 logger.info("running for accuracy %s%%" % (self.err_goal*100)) 884 self.gen_events = False 885 elif self.err_goal >= 1: 886 logger.info("Generating %s unweigthed events." % self.err_goal) 887 self.gen_events = True 888 self.err_goal = self.err_goal * self.gen_events_security # security
889
890 - def handle_seed(self):
891 """not needed but for gridpack --which is not handle here for the moment""" 892 return
893 894
895 - def find_job_for_event(self):
896 """return the list of channel that need to be improved""" 897 898 assert self.err_goal >=1 899 self.err_goal = int(self.err_goal) 900 901 goal_lum = self.err_goal/(self.results.axsec+1e-99) #pb^-1 902 logger.info('Effective Luminosity %s pb^-1', goal_lum) 903 904 all_channels = sum([list(P) for P in self.results],[]) 905 all_channels.sort(cmp= lambda x,y: 1 if y.get('luminosity') - \ 906 x.get('luminosity') > 0 else -1) 907 908 to_refine = [] 909 for C in all_channels: 910 if C.get('axsec') == 0: 911 continue 912 if goal_lum/(C.get('luminosity')+1e-99) >= 1 + (self.gen_events_security-1)/2: 913 logger.debug("channel %s is at %s (%s) (%s pb)", C.name, C.get('luminosity'), goal_lum/(C.get('luminosity')+1e-99), C.get('xsec')) 914 to_refine.append(C) 915 elif C.get('xerr') > max(C.get('axsec'), 916 (1/(100*math.sqrt(self.err_goal)))*all_channels[-1].get('axsec')): 917 to_refine.append(C) 918 919 logger.info('need to improve %s channels' % len(to_refine)) 920 return goal_lum, to_refine
921
922 - def update_html(self):
923 """update the html from this object since it contains all the information""" 924 925 926 run = self.cmd.results.current['run_name'] 927 if not os.path.exists(pjoin(self.cmd.me_dir, 'HTML', run)): 928 os.mkdir(pjoin(self.cmd.me_dir, 'HTML', run)) 929 930 unit = self.cmd.results.unit 931 P_text = "" 932 if self.results: 933 Presults = self.results 934 else: 935 self.results = sum_html.collect_result(self.cmd, None) 936 Presults = self.results 937 938 for P_comb in Presults: 939 P_text += P_comb.get_html(run, unit, self.cmd.me_dir) 940 941 Presults.write_results_dat(pjoin(self.cmd.me_dir,'SubProcesses', 'results.dat')) 942 943 fsock = open(pjoin(self.cmd.me_dir, 'HTML', run, 'results.html'),'w') 944 fsock.write(sum_html.results_header) 945 fsock.write('%s <dl>' % Presults.get_html(run, unit, self.cmd.me_dir)) 946 fsock.write('%s </dl></body>' % P_text) 947 948 self.cmd.results.add_detail('cross', Presults.xsec) 949 self.cmd.results.add_detail('error', Presults.xerru) 950 951 return Presults.xsec, Presults.xerru
952
953 954 -class gen_ximprove_v4(gen_ximprove):
955 956 # some hardcoded value which impact the generation 957 gen_events_security = 1.2 # multiply the number of requested event by this number for security 958 combining_job = 0 # allow to run multiple channel in sequence 959 max_request_event = 1000 # split jobs if a channel if it needs more than that 960 max_event_in_iter = 5000 961 min_event_in_iter = 1000 962 max_splitting = 130 # maximum duplication of a given channel 963 min_iter = 3 964 max_iter = 9 965 keep_grid_for_refine = False # only apply if needed to split the job 966 967 968
969 - def __init__(self, cmd, opt=None):
970 971 super(gen_ximprove_v4, self).__init__(cmd, opt) 972 973 if cmd.opts['accuracy'] < cmd._survey_options['accuracy'][1]: 974 self.increase_precision()
975
976 - def reset_multijob(self):
977 978 for path in misc.glob(pjoin('*', '*','multijob.dat'), pjoin(self.me_dir, 'SubProcesses')): 979 open(path,'w').write('0\n')
980
981 - def write_multijob(self, Channel, nb_split):
982 """ """ 983 if nb_split <=1: 984 return 985 f = open(pjoin(self.me_dir, 'SubProcesses', Channel.get('name'), 'multijob.dat'), 'w') 986 f.write('%i\n' % nb_split) 987 f.close()
988
989 - def increase_precision(self):
990 991 self.max_event_in_iter = 20000 992 self.min_events = 7500 993 if int(self.nhel) == 1: 994 self.min_event_in_iter *= 2**(self.cmd.proc_characteristics['nexternal']//3) 995 self.max_event_in_iter *= 2**(self.cmd.proc_characteristics['nexternal']//2) 996 997 self.gen_events_security = 1.3
998 999 alphabet = "abcdefghijklmnopqrstuvwxyz"
1000 - def get_job_for_event(self):
1001 """generate the script in order to generate a given number of event""" 1002 # correspond to write_gen in the fortran version 1003 1004 1005 goal_lum, to_refine = self.find_job_for_event() 1006 1007 #reset the potential multijob of previous run 1008 self.reset_multijob() 1009 1010 jobs = [] # list of the refine if some job are split is list of 1011 # dict with the parameter of the run. 1012 1013 # try to have a smart load on the cluster (not really important actually) 1014 if self.combining_job >1: 1015 # add a nice ordering for the jobs 1016 new_order = [] 1017 if self.combining_job % 2 == 0: 1018 for i in range(len(to_refine) //2): 1019 new_order.append(to_refine[i]) 1020 new_order.append(to_refine[-i-1]) 1021 if len(to_refine) % 2: 1022 new_order.append(to_refine[i+1]) 1023 else: 1024 for i in range(len(to_refine) //3): 1025 new_order.append(to_refine[i]) 1026 new_order.append(to_refine[-2*i-1]) 1027 new_order.append(to_refine[-2*i-2]) 1028 if len(to_refine) % 3 == 1: 1029 new_order.append(to_refine[i+1]) 1030 elif len(to_refine) % 3 == 2: 1031 new_order.append(to_refine[i+2]) 1032 #ensure that the reordering is done nicely 1033 assert set([id(C) for C in to_refine]) == set([id(C) for C in new_order]) 1034 to_refine = new_order 1035 1036 1037 # loop over the channel to refine 1038 for C in to_refine: 1039 #1. Compute the number of points are needed to reach target 1040 needed_event = goal_lum*C.get('axsec') 1041 nb_split = int(max(1,((needed_event-1)// self.max_request_event) +1)) 1042 if not self.split_channels: 1043 nb_split = 1 1044 if nb_split > self.max_splitting: 1045 nb_split = self.max_splitting 1046 nb_split=max(1, nb_split) 1047 1048 1049 #2. estimate how many points we need in each iteration 1050 if C.get('nunwgt') > 0: 1051 nevents = needed_event / nb_split * (C.get('nevents') / C.get('nunwgt')) 1052 #split by iter 1053 nevents = int(nevents / (2**self.min_iter-1)) 1054 else: 1055 nevents = self.max_event_in_iter 1056 1057 if nevents < self.min_event_in_iter: 1058 nb_split = int(nb_split * nevents / self.min_event_in_iter) + 1 1059 nevents = self.min_event_in_iter 1060 # 1061 # forbid too low/too large value 1062 nevents = max(self.min_event_in_iter, min(self.max_event_in_iter, nevents)) 1063 logger.debug("%s : need %s event. Need %s split job of %s points", C.name, needed_event, nb_split, nevents) 1064 1065 1066 # write the multi-job information 1067 self.write_multijob(C, nb_split) 1068 1069 packet = cluster.Packet((C.parent_name, C.name), 1070 combine_runs.CombineRuns, 1071 (pjoin(self.me_dir, 'SubProcesses', C.parent_name)), 1072 {"subproc": C.name, "nb_split":nb_split}) 1073 1074 1075 #create the info dict assume no splitting for the default 1076 info = {'name': self.cmd.results.current['run_name'], 1077 'script_name': 'unknown', 1078 'directory': C.name, # need to be change for splitted job 1079 'P_dir': C.parent_name, 1080 'Ppath': pjoin(self.cmd.me_dir, 'SubProcesses', C.parent_name), 1081 'offset': 1, # need to be change for splitted job 1082 'nevents': nevents, 1083 'maxiter': self.max_iter, 1084 'miniter': self.min_iter, 1085 'precision': -goal_lum/nb_split, 1086 'nhel': self.run_card['nhel'], 1087 'channel': C.name.replace('G',''), 1088 'grid_refinment' : 0, #no refinment of the grid 1089 'base_directory': '', #should be change in splitted job if want to keep the grid 1090 'packet': packet, 1091 } 1092 1093 if nb_split == 1: 1094 jobs.append(info) 1095 else: 1096 for i in range(nb_split): 1097 new_info = dict(info) 1098 new_info['offset'] = i+1 1099 new_info['directory'] += self.alphabet[i % 26] + str((i+1)//26) 1100 if self.keep_grid_for_refine: 1101 new_info['base_directory'] = info['directory'] 1102 jobs.append(new_info) 1103 1104 self.create_ajob(pjoin(self.me_dir, 'SubProcesses', 'refine.sh'), jobs)
1105 1106
1107 - def create_ajob(self, template, jobs, write_dir=None):
1108 """create the ajob""" 1109 1110 if not jobs: 1111 return 1112 1113 if not write_dir: 1114 write_dir = pjoin(self.me_dir, 'SubProcesses') 1115 1116 #filter the job according to their SubProcess directory # no mix submition 1117 P2job= collections.defaultdict(list) 1118 for j in jobs: 1119 P2job[j['P_dir']].append(j) 1120 if len(P2job) >1: 1121 for P in P2job.values(): 1122 self.create_ajob(template, P, write_dir) 1123 return 1124 1125 1126 #Here we can assume that all job are for the same directory. 1127 path = pjoin(write_dir, jobs[0]['P_dir']) 1128 1129 template_text = open(template, 'r').read() 1130 # special treatment if needed to combine the script 1131 # computes how many submition miss one job 1132 if self.combining_job > 1: 1133 skip1=0 1134 n_channels = len(jobs) 1135 nb_sub = n_channels // self.combining_job 1136 nb_job_in_last = n_channels % self.combining_job 1137 if nb_sub == 0: 1138 nb_sub = 1 1139 nb_job_in_last =0 1140 if nb_job_in_last: 1141 nb_sub +=1 1142 skip1 = self.combining_job - nb_job_in_last 1143 if skip1 > nb_sub: 1144 self.combining_job -=1 1145 return self.create_ajob(template, jobs, write_dir) 1146 combining_job = self.combining_job 1147 else: 1148 #define the variable for combining jobs even in not combine mode 1149 #such that we can use the same routine 1150 skip1=0 1151 combining_job =1 1152 nb_sub = len(jobs) 1153 1154 1155 nb_use = 0 1156 for i in range(nb_sub): 1157 script_number = i+1 1158 if i < skip1: 1159 nb_job = combining_job -1 1160 else: 1161 nb_job = min(combining_job, len(jobs)) 1162 fsock = open(pjoin(path, 'ajob%i' % script_number), 'w') 1163 for j in range(nb_use, nb_use + nb_job): 1164 if j> len(jobs): 1165 break 1166 info = jobs[j] 1167 info['script_name'] = 'ajob%i' % script_number 1168 info['keeplog'] = 'false' 1169 if "base_directory" not in info: 1170 info["base_directory"] = "./" 1171 fsock.write(template_text % info) 1172 nb_use += nb_job 1173 1174 fsock.close() 1175 return script_number
1176
1177 - def get_job_for_precision(self):
1178 """create the ajob to achieve a give precision on the total cross-section""" 1179 1180 1181 assert self.err_goal <=1 1182 xtot = abs(self.results.xsec) 1183 logger.info("Working on precision: %s %%" %(100*self.err_goal)) 1184 all_channels = sum([list(P) for P in self.results if P.mfactor],[]) 1185 limit = self.err_goal * xtot / len(all_channels) 1186 to_refine = [] 1187 rerr = 0 #error of the job not directly selected 1188 for C in all_channels: 1189 cerr = C.mfactor*(C.xerru + len(all_channels)*C.xerrc) 1190 if cerr > abs(limit): 1191 to_refine.append(C) 1192 else: 1193 rerr += cerr 1194 rerr *=rerr 1195 if not len(to_refine): 1196 return 1197 1198 # change limit since most don't contribute 1199 limit = math.sqrt((self.err_goal * xtot)**2 - rerr/math.sqrt(len(to_refine))) 1200 for C in to_refine[:]: 1201 cerr = C.mfactor*(C.xerru + len(to_refine)*C.xerrc) 1202 if cerr < limit: 1203 to_refine.remove(C) 1204 1205 # all the channel are now selected. create the channel information 1206 logger.info('need to improve %s channels' % len(to_refine)) 1207 1208 1209 jobs = [] # list of the refine if some job are split is list of 1210 # dict with the parameter of the run. 1211 1212 # loop over the channel to refine 1213 for C in to_refine: 1214 1215 #1. Determine how many events we need in each iteration 1216 yerr = C.mfactor*(C.xerru+len(to_refine)*C.xerrc) 1217 nevents = 0.2*C.nevents*(yerr/limit)**2 1218 1219 nb_split = int((nevents*(C.nunwgt/C.nevents)/self.max_request_event/ (2**self.min_iter-1))**(2/3)) 1220 nb_split = max(nb_split, 1) 1221 # **(2/3) to slow down the increase in number of jobs 1222 if nb_split > self.max_splitting: 1223 nb_split = self.max_splitting 1224 1225 if nb_split >1: 1226 nevents = nevents / nb_split 1227 self.write_multijob(C, nb_split) 1228 # forbid too low/too large value 1229 nevents = min(self.min_event_in_iter, max(self.max_event_in_iter, nevents)) 1230 1231 1232 #create the info dict assume no splitting for the default 1233 info = {'name': self.cmd.results.current['run_name'], 1234 'script_name': 'unknown', 1235 'directory': C.name, # need to be change for splitted job 1236 'P_dir': C.parent_name, 1237 'Ppath': pjoin(self.cmd.me_dir, 'SubProcesses', C.parent_name), 1238 'offset': 1, # need to be change for splitted job 1239 'nevents': nevents, 1240 'maxiter': self.max_iter, 1241 'miniter': self.min_iter, 1242 'precision': yerr/math.sqrt(nb_split)/(C.get('xsec')+ yerr), 1243 'nhel': self.run_card['nhel'], 1244 'channel': C.name.replace('G',''), 1245 'grid_refinment' : 1 1246 } 1247 1248 if nb_split == 1: 1249 jobs.append(info) 1250 else: 1251 for i in range(nb_split): 1252 new_info = dict(info) 1253 new_info['offset'] = i+1 1254 new_info['directory'] += self.alphabet[i % 26] + str((i+1)//26) 1255 jobs.append(new_info) 1256 self.create_ajob(pjoin(self.me_dir, 'SubProcesses', 'refine.sh'), jobs)
1257
1258 - def update_html(self):
1259 """update the html from this object since it contains all the information""" 1260 1261 1262 run = self.cmd.results.current['run_name'] 1263 if not os.path.exists(pjoin(self.cmd.me_dir, 'HTML', run)): 1264 os.mkdir(pjoin(self.cmd.me_dir, 'HTML', run)) 1265 1266 unit = self.cmd.results.unit 1267 P_text = "" 1268 if self.results: 1269 Presults = self.results 1270 else: 1271 self.results = sum_html.collect_result(self.cmd, None) 1272 Presults = self.results 1273 1274 for P_comb in Presults: 1275 P_text += P_comb.get_html(run, unit, self.cmd.me_dir) 1276 1277 Presults.write_results_dat(pjoin(self.cmd.me_dir,'SubProcesses', 'results.dat')) 1278 1279 fsock = open(pjoin(self.cmd.me_dir, 'HTML', run, 'results.html'),'w') 1280 fsock.write(sum_html.results_header) 1281 fsock.write('%s <dl>' % Presults.get_html(run, unit, self.cmd.me_dir)) 1282 fsock.write('%s </dl></body>' % P_text) 1283 1284 self.cmd.results.add_detail('cross', Presults.xsec) 1285 self.cmd.results.add_detail('error', Presults.xerru) 1286 1287 return Presults.xsec, Presults.xerru
1288
1289 1290 1291 1292 -class gen_ximprove_v4_nogridupdate(gen_ximprove_v4):
1293 1294 # some hardcoded value which impact the generation 1295 gen_events_security = 1.1 # multiply the number of requested event by this number for security 1296 combining_job = 0 # allow to run multiple channel in sequence 1297 max_request_event = 400 # split jobs if a channel if it needs more than that 1298 max_event_in_iter = 500 1299 min_event_in_iter = 250 1300 max_splitting = 260 # maximum duplication of a given channel 1301 min_iter = 2 1302 max_iter = 6 1303 keep_grid_for_refine = True 1304 1305
1306 - def __init__(self, cmd, opt=None):
1307 1308 gen_ximprove.__init__(cmd, opt) 1309 1310 if cmd.proc_characteristics['loopinduced'] and \ 1311 cmd.proc_characteristics['nexternal'] > 2: 1312 self.increase_parralelization(cmd.proc_characteristics['nexternal'])
1313
1314 - def increase_parralelization(self, nexternal):
1315 1316 self.max_splitting = 1000 1317 1318 if self.run_card['refine_evt_by_job'] != -1: 1319 pass 1320 elif nexternal == 3: 1321 self.max_request_event = 200 1322 elif nexternal == 4: 1323 self.max_request_event = 100 1324 elif nexternal >= 5: 1325 self.max_request_event = 50 1326 self.min_event_in_iter = 125 1327 self.max_iter = 5
1328
1329 -class gen_ximprove_share(gen_ximprove, gensym):
1330 """Doing the refine in multicore. Each core handle a couple of PS point.""" 1331 1332 nb_ps_by_job = 2000 1333 mode = "refine" 1334 gen_events_security = 1.15 1335 # Note the real security is lower since we stop the jobs if they are at 96% 1336 # of this target. 1337
1338 - def __init__(self, *args, **opts):
1339 1340 super(gen_ximprove_share, self).__init__(*args, **opts) 1341 self.generated_events = {} 1342 self.splitted_for_dir = lambda x,y : self.splitted_Pdir[(x,y)]
1343 1344
1345 - def get_job_for_event(self):
1346 """generate the script in order to generate a given number of event""" 1347 # correspond to write_gen in the fortran version 1348 1349 1350 goal_lum, to_refine = self.find_job_for_event() 1351 self.goal_lum = goal_lum 1352 1353 # loop over the channel to refine to find the number of PS point to launch 1354 total_ps_points = 0 1355 channel_to_ps_point = [] 1356 for C in to_refine: 1357 #0. remove previous events files 1358 try: 1359 os.remove(pjoin(self.me_dir, "SubProcesses",C.parent_name, C.name, "events.lhe")) 1360 except: 1361 pass 1362 1363 #1. Compute the number of points are needed to reach target 1364 needed_event = goal_lum*C.get('axsec') 1365 if needed_event == 0: 1366 continue 1367 #2. estimate how many points we need in each iteration 1368 if C.get('nunwgt') > 0: 1369 nevents = needed_event * (C.get('nevents') / C.get('nunwgt')) 1370 #split by iter 1371 nevents = int(nevents / (2**self.min_iter-1)) 1372 else: 1373 nb_split = int(max(1,((needed_event-1)// self.max_request_event) +1)) 1374 if not self.split_channels: 1375 nb_split = 1 1376 if nb_split > self.max_splitting: 1377 nb_split = self.max_splitting 1378 nevents = self.max_event_in_iter * self.max_splitting 1379 else: 1380 nevents = self.max_event_in_iter * nb_split 1381 1382 if nevents > self.max_splitting*self.max_event_in_iter: 1383 logger.warning("Channel %s/%s has a very low efficiency of unweighting. Might not be possible to reach target" % \ 1384 (C.name, C.parent_name)) 1385 nevents = self.max_event_in_iter * self.max_splitting 1386 1387 total_ps_points += nevents 1388 channel_to_ps_point.append((C, nevents)) 1389 1390 if self.cmd.options["run_mode"] == 1: 1391 if self.cmd.options["cluster_size"]: 1392 nb_ps_by_job = total_ps_points /int(self.cmd.options["cluster_size"]) 1393 else: 1394 nb_ps_by_job = self.nb_ps_by_job 1395 elif self.cmd.options["run_mode"] == 2: 1396 remain = total_ps_points % self.cmd.options["nb_core"] 1397 if remain: 1398 nb_ps_by_job = 1 + (total_ps_points - remain) / self.cmd.options["nb_core"] 1399 else: 1400 nb_ps_by_job = total_ps_points / self.cmd.options["nb_core"] 1401 else: 1402 nb_ps_by_job = self.nb_ps_by_job 1403 1404 nb_ps_by_job = int(max(nb_ps_by_job, 500)) 1405 1406 for C, nevents in channel_to_ps_point: 1407 if nevents % nb_ps_by_job: 1408 nb_job = 1 + int(nevents // nb_ps_by_job) 1409 else: 1410 nb_job = int(nevents // nb_ps_by_job) 1411 submit_ps = min(nevents, nb_ps_by_job) 1412 if nb_job == 1: 1413 submit_ps = max(submit_ps, self.min_event_in_iter) 1414 self.create_resubmit_one_iter(C.parent_name, C.name[1:], submit_ps, nb_job, step=0) 1415 needed_event = goal_lum*C.get('xsec') 1416 logger.debug("%s/%s : need %s event. Need %s split job of %s points", C.parent_name, C.name, needed_event, nb_job, submit_ps)
1417 1418
1419 - def combine_iteration(self, Pdir, G, step):
1420 1421 grid_calculator, cross, error = self.combine_grid(Pdir, G, step) 1422 1423 # collect all the generated_event 1424 Gdirs = [] #build the the list of directory 1425 for i in range(self.splitted_for_dir(Pdir, G)): 1426 path = pjoin(Pdir, "G%s_%s" % (G, i+1)) 1427 Gdirs.append(path) 1428 assert len(grid_calculator.results) == len(Gdirs) == self.splitted_for_dir(Pdir, G) 1429 1430 1431 # Check how many events are going to be kept after un-weighting. 1432 needed_event = cross * self.goal_lum 1433 if needed_event == 0: 1434 return 0 1435 # check that the number of events requested is not higher than the actual 1436 # total number of events to generate. 1437 if self.err_goal >=1: 1438 if needed_event > self.gen_events_security * self.err_goal: 1439 needed_event = int(self.gen_events_security * self.err_goal) 1440 1441 if (Pdir, G) in self.generated_events: 1442 old_nunwgt, old_maxwgt = self.generated_events[(Pdir, G)] 1443 else: 1444 old_nunwgt, old_maxwgt = 0, 0 1445 1446 if old_nunwgt == 0 and os.path.exists(pjoin(Pdir,"G%s" % G, "events.lhe")): 1447 # possible for second refine. 1448 lhe = lhe_parser.EventFile(pjoin(Pdir,"G%s" % G, "events.lhe")) 1449 old_nunwgt = lhe.unweight(None, trunc_error=0.005, log_level=0) 1450 old_maxwgt = lhe.max_wgt 1451 1452 1453 1454 maxwgt = max(grid_calculator.get_max_wgt(), old_maxwgt) 1455 new_evt = grid_calculator.get_nunwgt(maxwgt) 1456 efficiency = new_evt / sum([R.nevents for R in grid_calculator.results]) 1457 nunwgt = old_nunwgt * old_maxwgt / maxwgt 1458 nunwgt += new_evt 1459 1460 # check the number of event for this iteration alone 1461 one_iter_nb_event = max(grid_calculator.get_nunwgt(),1) 1462 drop_previous_iteration = False 1463 # compare the number of events to generate if we discard the previous iteration 1464 n_target_one_iter = (needed_event-one_iter_nb_event) / ( one_iter_nb_event/ sum([R.nevents for R in grid_calculator.results])) 1465 n_target_combined = (needed_event-nunwgt) / efficiency 1466 if n_target_one_iter < n_target_combined: 1467 # the last iteration alone has more event that the combine iteration. 1468 # it is therefore interesting to drop previous iteration. 1469 drop_previous_iteration = True 1470 nunwgt = one_iter_nb_event 1471 maxwgt = grid_calculator.get_max_wgt() 1472 new_evt = nunwgt 1473 efficiency = ( one_iter_nb_event/ sum([R.nevents for R in grid_calculator.results])) 1474 1475 try: 1476 if drop_previous_iteration: 1477 raise IOError 1478 output_file = open(pjoin(Pdir,"G%s" % G, "events.lhe"), 'a') 1479 except IOError: 1480 output_file = open(pjoin(Pdir,"G%s" % G, "events.lhe"), 'w') 1481 1482 misc.call(["cat"] + [pjoin(d, "events.lhe") for d in Gdirs], 1483 stdout=output_file) 1484 output_file.close() 1485 # For large number of iteration. check the number of event by doing the 1486 # real unweighting. 1487 if nunwgt < 0.6 * needed_event and step > self.min_iter: 1488 lhe = lhe_parser.EventFile(output_file.name) 1489 old_nunwgt =nunwgt 1490 nunwgt = lhe.unweight(None, trunc_error=0.01, log_level=0) 1491 1492 1493 self.generated_events[(Pdir, G)] = (nunwgt, maxwgt) 1494 1495 # misc.sprint("Adding %s event to %s. Currently at %s" % (new_evt, G, nunwgt)) 1496 # check what to do 1497 if nunwgt >= int(0.96*needed_event)+1: # 0.96*1.15=1.10 =real security 1498 # We did it. 1499 logger.info("found enough event for %s/G%s" % (os.path.basename(Pdir), G)) 1500 self.write_results(grid_calculator, cross, error, Pdir, G, step, efficiency) 1501 return 0 1502 elif step >= self.max_iter: 1503 logger.debug("fail to find enough event") 1504 self.write_results(grid_calculator, cross, error, Pdir, G, step, efficiency) 1505 return 0 1506 1507 nb_split_before = len(grid_calculator.results) 1508 nevents = grid_calculator.results[0].nevents 1509 if nevents == 0: # possible if some integral returns 0 1510 nevents = max(g.nevents for g in grid_calculator.results) 1511 1512 need_ps_point = (needed_event - nunwgt)/(efficiency+1e-99) 1513 need_job = need_ps_point // nevents + 1 1514 1515 if step < self.min_iter: 1516 # This is normal but check if we are on the good track 1517 job_at_first_iter = nb_split_before/2**(step-1) 1518 expected_total_job = job_at_first_iter * (2**self.min_iter-1) 1519 done_job = job_at_first_iter * (2**step-1) 1520 expected_remaining_job = expected_total_job - done_job 1521 1522 logger.debug("efficiency status (smaller is better): %s", need_job/expected_remaining_job) 1523 # increase if needed but not too much 1524 need_job = min(need_job, expected_remaining_job*1.25) 1525 1526 nb_job = (need_job-0.5)//(2**(self.min_iter-step)-1) + 1 1527 nb_job = max(1, nb_job) 1528 grid_calculator.write_grid_for_submission(Pdir,G, 1529 self.splitted_for_dir(Pdir, G), nb_job*nevents ,mode=self.mode, 1530 conservative_factor=self.max_iter) 1531 logger.info("%s/G%s is at %i/%i (%.2g%%) event. Resubmit %i job at iteration %i." \ 1532 % (os.path.basename(Pdir), G, int(nunwgt),int(needed_event)+1, 1533 (float(nunwgt)/needed_event)*100.0 if needed_event>0.0 else 0.0, 1534 nb_job, step)) 1535 self.create_resubmit_one_iter(Pdir, G, nevents, nb_job, step) 1536 #self.create_job(Pdir, G, nb_job, nevents, step) 1537 1538 elif step < self.max_iter: 1539 if step + 1 == self.max_iter: 1540 need_job = 1.20 * need_job # avoid to have just too few event. 1541 1542 nb_job = int(min(need_job, nb_split_before*1.5)) 1543 grid_calculator.write_grid_for_submission(Pdir,G, 1544 self.splitted_for_dir(Pdir, G), nb_job*nevents ,mode=self.mode, 1545 conservative_factor=self.max_iter) 1546 1547 1548 logger.info("%s/G%s is at %i/%i ('%.2g%%') event. Resubmit %i job at iteration %i." \ 1549 % (os.path.basename(Pdir), G, int(nunwgt),int(needed_event)+1, 1550 (float(nunwgt)/needed_event)*100.0 if needed_event>0.0 else 0.0, 1551 nb_job, step)) 1552 self.create_resubmit_one_iter(Pdir, G, nevents, nb_job, step) 1553 1554 1555 1556 return 0
1557 1558
1559 - def write_results(self, grid_calculator, cross, error, Pdir, G, step, efficiency):
1560 1561 #compute the value 1562 if cross == 0: 1563 abscross,nw, luminosity = 0, 0, 0 1564 wgt, maxit,nunwgt, wgt, nevents = 0,0,0,0,0 1565 error = 0 1566 else: 1567 grid_calculator.results.compute_values() 1568 abscross = self.abscross[(Pdir,G)]/self.sigma[(Pdir,G)] 1569 nunwgt, wgt = self.generated_events[(Pdir, G)] 1570 nw = int(nunwgt / efficiency) 1571 nunwgt = int(nunwgt) 1572 maxit = step 1573 nevents = nunwgt 1574 # make the unweighting to compute the number of events: 1575 luminosity = nunwgt/cross 1576 1577 #format the results.dat 1578 def fstr(nb): 1579 data = '%E' % nb 1580 nb, power = data.split('E') 1581 nb = float(nb) /10 1582 power = int(power) + 1 1583 return '%.5fE%+03i' %(nb,power)
1584 line = '%s %s %s %i %i %i %i %s %s %s 0.0 0.0 0\n' % \ 1585 (fstr(cross), fstr(error*cross), fstr(error*cross), 1586 nevents, nw, maxit,nunwgt, 1587 fstr(luminosity), fstr(wgt), fstr(abscross)) 1588 1589 fsock = open(pjoin(self.me_dir,'SubProcesses' , Pdir, 'G%s' % G, 1590 'results.dat'),'w') 1591 fsock.writelines(line) 1592 fsock.close()
1593
1594 1595 1596 1597 -class gen_ximprove_gridpack(gen_ximprove_v4):
1598 1599 min_iter = 1 1600 max_iter = 13 1601 max_request_event = 1e12 # split jobs if a channel if it needs more than that 1602 max_event_in_iter = 4000 1603 min_event_in_iter = 500 1604 combining_job = sys.maxint 1605 gen_events_security = 1.00 1606
1607 - def __new__(cls, *args, **opts):
1608 1609 return super(gen_ximprove, cls).__new__(cls, *args, **opts)
1610
1611 - def __init__(self, *args, **opts):
1612 1613 self.ngran = -1 1614 self.gscalefact = {} 1615 self.readonly = False 1616 if 'ngran' in opts: 1617 self.gran = opts['ngran'] 1618 # del opts['ngran'] 1619 if 'readonly' in opts: 1620 self.readonly = opts['readonly'] 1621 super(gen_ximprove_gridpack,self).__init__(*args, **opts) 1622 if self.ngran == -1: 1623 self.ngran = 1
1624
1625 - def find_job_for_event(self):
1626 """return the list of channel that need to be improved""" 1627 import random 1628 1629 assert self.err_goal >=1 1630 self.err_goal = int(self.err_goal) 1631 self.gscalefact = {} 1632 1633 xtot = self.results.axsec 1634 goal_lum = self.err_goal/(xtot+1e-99) #pb^-1 1635 # logger.info('Effective Luminosity %s pb^-1', goal_lum) 1636 1637 all_channels = sum([list(P) for P in self.results],[]) 1638 all_channels.sort(cmp= lambda x,y: 1 if y.get('luminosity') - \ 1639 x.get('luminosity') > 0 else -1) 1640 1641 to_refine = [] 1642 for C in all_channels: 1643 tag = C.get('name') 1644 self.gscalefact[tag] = 0 1645 R = random.random() 1646 if C.get('axsec') == 0: 1647 continue 1648 if (goal_lum * C.get('axsec') < R*self.ngran ): 1649 continue # no event to generate events 1650 self.gscalefact[tag] = max(1, 1/(goal_lum * C.get('axsec')/ self.ngran)) 1651 #need to generate events 1652 logger.debug('request events for ', C.get('name'), 'cross=', 1653 C.get('axsec'), 'needed events = ', goal_lum * C.get('axsec')) 1654 to_refine.append(C) 1655 1656 logger.info('need to improve %s channels' % len(to_refine)) 1657 return goal_lum, to_refine
1658
1659 - def get_job_for_event(self):
1660 """generate the script in order to generate a given number of event""" 1661 # correspond to write_gen in the fortran version 1662 1663 1664 goal_lum, to_refine = self.find_job_for_event() 1665 1666 jobs = [] # list of the refine if some job are split is list of 1667 # dict with the parameter of the run. 1668 1669 # loop over the channel to refine 1670 for C in to_refine: 1671 #1. Compute the number of points are needed to reach target 1672 needed_event = max(goal_lum*C.get('axsec'), self.ngran) 1673 nb_split = 1 1674 1675 #2. estimate how many points we need in each iteration 1676 if C.get('nunwgt') > 0: 1677 nevents = needed_event / nb_split * (C.get('nevents') / C.get('nunwgt')) 1678 #split by iter 1679 nevents = int(nevents / (2**self.min_iter-1)) 1680 else: 1681 nevents = self.max_event_in_iter 1682 1683 if nevents < self.min_event_in_iter: 1684 nevents = self.min_event_in_iter 1685 # 1686 # forbid too low/too large value 1687 nevents = max(self.min_event_in_iter, min(self.max_event_in_iter, nevents)) 1688 logger.debug("%s : need %s event. Need %s split job of %s points", C.name, needed_event, nb_split, nevents) 1689 1690 1691 #create the info dict assume no splitting for the default 1692 info = {'name': self.cmd.results.current['run_name'], 1693 'script_name': 'unknown', 1694 'directory': C.name, # need to be change for splitted job 1695 'P_dir': os.path.basename(C.parent_name), 1696 'offset': 1, # need to be change for splitted job 1697 'Ppath': pjoin(self.cmd.me_dir, 'SubProcesses', C.parent_name), 1698 'nevents': nevents, #int(nevents*self.gen_events_security)+1, 1699 'maxiter': self.max_iter, 1700 'miniter': self.min_iter, 1701 'precision': -1*int(needed_event)/C.get('axsec'), 1702 'requested_event': needed_event, 1703 'nhel': self.run_card['nhel'], 1704 'channel': C.name.replace('G',''), 1705 'grid_refinment' : 0, #no refinment of the grid 1706 'base_directory': '', #should be change in splitted job if want to keep the grid 1707 'packet': None, 1708 } 1709 1710 1711 jobs.append(info) 1712 1713 1714 write_dir = '.' if self.readonly else None 1715 self.create_ajob(pjoin(self.me_dir, 'SubProcesses', 'refine.sh'), jobs, write_dir) 1716 1717 done = [] 1718 for j in jobs: 1719 if j['P_dir'] in done: 1720 continue 1721 done.append(j['P_dir']) 1722 # set the working directory path. 1723 pwd = pjoin(os.getcwd(),j['P_dir']) if self.readonly else pjoin(self.me_dir, 'SubProcesses', j['P_dir']) 1724 exe = pjoin(pwd, 'ajob1') 1725 st = os.stat(exe) 1726 os.chmod(exe, st.st_mode | stat.S_IEXEC) 1727 1728 # run the code\ 1729 cluster.onecore.launch_and_wait(exe, cwd=pwd, packet_member=j['packet']) 1730 write_dir = '.' if self.readonly else pjoin(self.me_dir, 'SubProcesses') 1731 1732 self.check_events(goal_lum, to_refine, jobs, write_dir)
1733
1734 - def check_events(self, goal_lum, to_refine, jobs, Sdir):
1735 """check that we get the number of requested events if not resubmit.""" 1736 1737 new_jobs = [] 1738 1739 for C, job_info in zip(to_refine, jobs): 1740 P = job_info['P_dir'] 1741 G = job_info['channel'] 1742 axsec = C.get('axsec') 1743 requested_events= job_info['requested_event'] 1744 1745 1746 new_results = sum_html.OneResult((P,G)) 1747 new_results.read_results(pjoin(Sdir,P, 'G%s'%G, 'results.dat')) 1748 1749 # need to resubmit? 1750 if new_results.get('nunwgt') < requested_events: 1751 pwd = pjoin(os.getcwd(),job_info['P_dir'],'G%s'%G) if self.readonly else \ 1752 pjoin(self.me_dir, 'SubProcesses', job_info['P_dir'],'G%s'%G) 1753 job_info['requested_event'] -= new_results.get('nunwgt') 1754 job_info['precision'] -= -1*job_info['requested_event']/axsec 1755 job_info['offset'] += 1 1756 new_jobs.append(job_info) 1757 files.mv(pjoin(pwd, 'events.lhe'), pjoin(pwd, 'events.lhe.previous')) 1758 1759 if new_jobs: 1760 self.create_ajob(pjoin(self.me_dir, 'SubProcesses', 'refine.sh'), new_jobs, Sdir) 1761 1762 done = [] 1763 for j in new_jobs: 1764 if j['P_dir'] in done: 1765 continue 1766 G = j['channel'] 1767 # set the working directory path. 1768 pwd = pjoin(os.getcwd(),j['P_dir']) if self.readonly \ 1769 else pjoin(self.me_dir, 'SubProcesses', j['P_dir']) 1770 exe = pjoin(pwd, 'ajob1') 1771 st = os.stat(exe) 1772 os.chmod(exe, st.st_mode | stat.S_IEXEC) 1773 1774 # run the code 1775 cluster.onecore.launch_and_wait(exe, cwd=pwd, packet_member=j['packet']) 1776 pwd = pjoin(pwd, 'G%s'%G) 1777 # concatanate with old events file 1778 files.put_at_end(pjoin(pwd, 'events.lhe'),pjoin(pwd, 'events.lhe.previous')) 1779 1780 return self.check_events(goal_lum, to_refine, new_jobs, Sdir)
1781