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