Package madgraph :: Package loop :: Module loop_exporters
[hide private]
[frames] | no frames]

Source Code for Module madgraph.loop.loop_exporters

   1  ################################################################################ 
   2  # 
   3  # Copyright (c) 2009 The MadGraph5_aMC@NLO Development team and Contributors 
   4  # 
   5  # This file is a part of the MadGraph5_aMC@NLO project, an application which  
   6  # automatically generates Feynman diagrams and matrix elements for arbitrary 
   7  # high-energy processes in the Standard Model and beyond. 
   8  # 
   9  # It is subject to the MadGraph5_aMC@NLO license which should accompany this  
  10  # distribution. 
  11  # 
  12  # For more information, visit madgraph.phys.ucl.ac.be and amcatnlo.web.cern.ch 
  13  # 
  14  ################################################################################ 
  15  """Methods and classes to export matrix elements to v4 format.""" 
  16   
  17  import copy 
  18  import fractions 
  19  import glob 
  20  import logging 
  21  import os 
  22  import sys 
  23  import re 
  24  import shutil 
  25  import subprocess 
  26  import itertools 
  27  import time 
  28  import datetime 
  29   
  30  import aloha 
  31   
  32  import madgraph.core.base_objects as base_objects 
  33  import madgraph.core.color_algebra as color 
  34  import madgraph.core.helas_objects as helas_objects 
  35  import madgraph.iolibs.drawing_eps as draw 
  36  import madgraph.iolibs.files as files 
  37  import madgraph.iolibs.group_subprocs as group_subprocs 
  38  import madgraph.various.misc as misc 
  39  import madgraph.various.q_polynomial as q_polynomial 
  40  import madgraph.iolibs.file_writers as writers 
  41  import madgraph.iolibs.gen_infohtml as gen_infohtml 
  42  import madgraph.iolibs.template_files as template_files 
  43  import madgraph.iolibs.ufo_expression_parsers as parsers 
  44  import madgraph.iolibs.export_v4 as export_v4 
  45  import madgraph.various.diagram_symmetry as diagram_symmetry 
  46  import madgraph.various.process_checks as process_checks 
  47  import madgraph.various.progressbar as pbar 
  48  import madgraph.core.color_amp as color_amp 
  49  import madgraph.iolibs.helas_call_writers as helas_call_writers 
  50  import models.check_param_card as check_param_card 
  51  from madgraph.loop.loop_base_objects import LoopDiagram 
  52  from madgraph.loop.MadLoopBannerStyles import MadLoopBannerStyles 
  53   
  54  import madgraph.various.banner as banner_mod 
  55   
  56  pjoin = os.path.join 
  57   
  58  import aloha.create_aloha as create_aloha 
  59  import models.write_param_card as param_writer 
  60  from madgraph import MadGraph5Error, MG5DIR, InvalidCmd 
  61  from madgraph.iolibs.files import cp, ln, mv 
  62  pjoin = os.path.join 
  63  _file_path = os.path.split(os.path.dirname(os.path.realpath(__file__)))[0] + '/' 
  64  logger = logging.getLogger('madgraph.loop_exporter') 
  65   
  66  #=============================================================================== 
  67  # LoopExporterFortran 
  68  #=============================================================================== 
69 -class LoopExporterFortran(object):
70 """ Class to define general helper functions to the different 71 loop fortran exporters (ME, SA, MEGroup, etc..) which will inherit both 72 from this class AND from the corresponding ProcessExporterFortran(ME,SA,...). 73 It plays the same role as ProcessExporterFrotran and simply defines here 74 loop-specific helpers functions necessary for all loop exporters. 75 Notice that we do not have LoopExporterFortran inheriting from 76 ProcessExporterFortran but give access to arguments like dir_path and 77 clean using options. This avoids method resolution object ambiguity""" 78 79 default_opt = {'clean': False, 'complex_mass':False, 80 'export_format':'madloop', 'mp':True, 81 'loop_dir':'', 'cuttools_dir':'', 82 'fortran_compiler':'gfortran', 83 'SubProc_prefix': 'P', 84 'output_dependencies': 'external', 85 'compute_color_flows': False, 86 'mode':''} 87 88
89 - def __init__(self, mgme_dir="", dir_path = "", opt=None):
90 """Initiate the LoopExporterFortran with directory information on where 91 to find all the loop-related source files, like CutTools""" 92 93 94 self.opt = dict(self.default_opt) 95 if opt: 96 self.opt.update(opt) 97 98 self.SubProc_prefix = self.opt['SubProc_prefix'] 99 self.loop_dir = self.opt['loop_dir'] 100 self.cuttools_dir = self.opt['cuttools_dir'] 101 self.fortran_compiler = self.opt['fortran_compiler'] 102 self.dependencies = self.opt['output_dependencies'] 103 self.compute_color_flows = self.opt['compute_color_flows'] 104 105 super(LoopExporterFortran,self).__init__(mgme_dir, dir_path, self.opt)
106 107 169
170 - def get_aloha_model(self, model):
171 """ Caches the aloha model created here as an attribute of the loop 172 exporter so that it can later be used in the LoopHelasMatrixElement 173 in the function compute_all_analytic_information for recycling aloha 174 computations across different LoopHelasMatrixElements steered by the 175 same loop exporter. 176 """ 177 if not hasattr(self, 'aloha_model'): 178 self.aloha_model = create_aloha.AbstractALOHAModel(model.get('name')) 179 return self.aloha_model
180 181 #=========================================================================== 182 # write the multiple-precision header files 183 #===========================================================================
184 - def write_mp_files(self, writer_mprec, writer_mpc):
185 """Write the cts_mprec.h and cts_mpc.h""" 186 187 file = open(os.path.join(self.cuttools_dir, 'src/cts/cts_mprec.h')).read() 188 writer_mprec.writelines(file) 189 190 file = open(os.path.join(self.cuttools_dir, 'src/cts/cts_mpc.h')).read() 191 file = file.replace('&','') 192 writer_mpc.writelines(file) 193 194 return True
195 196 #=============================================================================== 197 # LoopProcessExporterFortranSA 198 #===============================================================================
199 -class LoopProcessExporterFortranSA(LoopExporterFortran, 200 export_v4.ProcessExporterFortranSA):
201 202 """Class to take care of exporting a set of loop matrix elements in the 203 Fortran format.""" 204 205 template_dir=os.path.join(_file_path,'iolibs/template_files/loop') 206 madloop_makefile_name = 'makefile' 207 208 MadLoop_banner = MadLoopBannerStyles.get_MadLoop_Banner( 209 style='classic2', color='green', 210 top_frame_char = '=', bottom_frame_char = '=', 211 left_frame_char = '{',right_frame_char = '}', 212 print_frame=True, side_margin = 7, up_margin = 1) 213
214 - def copy_v4template(self, modelname = ''):
215 """Additional actions needed to setup the Template. 216 """ 217 super(LoopProcessExporterFortranSA, self).copy_v4template(modelname) 218 219 self.loop_additional_template_setup()
220
221 - def loop_additional_template_setup(self, copy_Source_makefile = True):
222 """ Perform additional actions specific for this class when setting 223 up the template with the copy_v4template function.""" 224 225 # We must change some files to their version for NLO computations 226 cpfiles= ["Cards/MadLoopParams.dat", 227 "SubProcesses/MadLoopParamReader.f", 228 "SubProcesses/MadLoopParams.inc"] 229 if copy_Source_makefile: 230 cpfiles.append("Source/makefile") 231 232 for file in cpfiles: 233 shutil.copy(os.path.join(self.loop_dir,'StandAlone/', file), 234 os.path.join(self.dir_path, file)) 235 236 self.MadLoopparam = banner_mod.MadLoopParam(pjoin(self.loop_dir,'StandAlone', 237 'Cards', 'MadLoopParams.dat')) 238 # write the output file 239 self.MadLoopparam.write(pjoin(self.dir_path,"SubProcesses", 240 "MadLoopParams.dat")) 241 242 # We might need to give a different name to the MadLoop makefile\ 243 shutil.copy(pjoin(self.loop_dir,'StandAlone','SubProcesses','makefile'), 244 pjoin(self.dir_path, 'SubProcesses',self.madloop_makefile_name)) 245 246 # Write SubProcesses/MadLoop_makefile_definitions with dummy variables 247 # for the non-optimized output 248 link_tir_libs=[] 249 tir_libs=[] 250 251 filePath = pjoin(self.dir_path, 'SubProcesses', 252 'MadLoop_makefile_definitions') 253 calls = self.write_loop_makefile_definitions( 254 writers.MakefileWriter(filePath),link_tir_libs,tir_libs) 255 256 257 # We need minimal editing of MadLoopCommons.f 258 MadLoopCommon = open(os.path.join(self.loop_dir,'StandAlone', 259 "SubProcesses","MadLoopCommons.inc")).read() 260 writer = writers.FortranWriter(os.path.join(self.dir_path, 261 "SubProcesses","MadLoopCommons.f")) 262 writer.writelines(MadLoopCommon%{ 263 'print_banner_commands':self.MadLoop_banner}) 264 writer.close() 265 266 267 # Copy the whole MadLoop5_resources directory (empty at this stage) 268 if not os.path.exists(pjoin(self.dir_path,'SubProcesses', 269 'MadLoop5_resources')): 270 cp(pjoin(self.loop_dir,'StandAlone','SubProcesses', 271 'MadLoop5_resources'),pjoin(self.dir_path,'SubProcesses')) 272 273 # Link relevant cards from Cards inside the MadLoop5_resources 274 ln(pjoin(self.dir_path,'SubProcesses','MadLoopParams.dat'), 275 pjoin(self.dir_path,'SubProcesses','MadLoop5_resources')) 276 ln(pjoin(self.dir_path,'Cards','param_card.dat'), 277 pjoin(self.dir_path,'SubProcesses','MadLoop5_resources')) 278 ln(pjoin(self.dir_path,'Cards','ident_card.dat'), 279 pjoin(self.dir_path,'SubProcesses','MadLoop5_resources')) 280 281 # And remove check_sa in the SubProcess folder since now there is a 282 # check_sa tailored to each subprocess. 283 if os.path.isfile(pjoin(self.dir_path,'SubProcesses','check_sa.f')): 284 os.remove(pjoin(self.dir_path,'SubProcesses','check_sa.f')) 285 286 cwd = os.getcwd() 287 dirpath = os.path.join(self.dir_path, 'SubProcesses') 288 try: 289 os.chdir(dirpath) 290 except os.error: 291 logger.error('Could not cd to directory %s' % dirpath) 292 return 0 293 294 # Write the cts_mpc.h and cts_mprec.h files imported from CutTools 295 self.write_mp_files(writers.FortranWriter('cts_mprec.h'),\ 296 writers.FortranWriter('cts_mpc.h')) 297 298 # Return to original PWD 299 os.chdir(cwd) 300 301 # We must link the CutTools to the Library folder of the active Template 302 super(LoopProcessExporterFortranSA, self).link_CutTools(self.dir_path)
303 304 # This function is placed here and not in optimized exporterd, 305 # because the same makefile.inc should be used in all cases.
306 - def write_loop_makefile_definitions(self, writer, link_tir_libs, 307 tir_libs,tir_include=[]):
308 """ Create the file makefile which links to the TIR libraries.""" 309 310 file = open(os.path.join(self.loop_dir,'StandAlone', 311 'SubProcesses','MadLoop_makefile_definitions.inc')).read() 312 replace_dict={} 313 replace_dict['link_tir_libs']=' '.join(link_tir_libs) 314 replace_dict['tir_libs']=' '.join(tir_libs) 315 replace_dict['dotf']='%.f' 316 replace_dict['prefix']= self.SubProc_prefix 317 replace_dict['doto']='%.o' 318 replace_dict['tir_include']=' '.join(tir_include) 319 file=file%replace_dict 320 if writer: 321 writer.writelines(file) 322 else: 323 return file
324
325 - def convert_model_to_mg4(self, model, wanted_lorentz = [], 326 wanted_couplings = []):
327 """ Caches the aloha model created here when writing out the aloha 328 fortran subroutine. 329 """ 330 self.get_aloha_model(model) 331 super(LoopProcessExporterFortranSA, self).convert_model_to_mg4(model, 332 wanted_lorentz = wanted_lorentz, wanted_couplings = wanted_couplings)
333
334 - def get_ME_identifier(self, matrix_element, 335 group_number = None, group_elem_number = None):
336 """ A function returning a string uniquely identifying the matrix 337 element given in argument so that it can be used as a prefix to all 338 MadLoop5 subroutines and common blocks related to it. This allows 339 to compile several processes into one library as requested by the 340 BLHA (Binoth LesHouches Accord) guidelines. 341 The arguments group_number and proc_id are just for the LoopInduced 342 output with MadEvent.""" 343 344 # When disabling the loop grouping in the LoopInduced MadEvent output, 345 # we have only the group_number set and the proc_id set to None. In this 346 # case we don't print the proc_id. 347 if (not group_number is None) and group_elem_number is None: 348 return 'ML5_%d_%s_'%(matrix_element.get('processes')[0].get('id'), 349 group_number) 350 elif group_number is None or group_elem_number is None: 351 return 'ML5_%d_'%matrix_element.get('processes')[0].get('id') 352 else: 353 return 'ML5_%d_%s_%s_'%(matrix_element.get('processes')[0].get('id'), 354 group_number, group_elem_number)
355
356 - def get_SubProc_folder_name(self, process, 357 group_number = None, group_elem_number = None):
358 """Returns the name of the SubProcess directory, which can contain 359 the process goup and group element number for the case of loop-induced 360 integration with MadEvent.""" 361 362 # When disabling the loop grouping in the LoopInduced MadEvent output, 363 # we have only the group_number set and the proc_id set to None. In this 364 # case we don't print the proc_id. 365 if not group_number is None and group_elem_number is None: 366 return "%s%d_%s_%s"%(self.SubProc_prefix, process.get('id'), 367 group_number,process.shell_string(print_id=False)) 368 elif group_number is None or group_elem_number is None: 369 return "%s%s" %(self.SubProc_prefix,process.shell_string()) 370 else: 371 return "%s%d_%s_%s_%s"%(self.SubProc_prefix, process.get('id'), 372 group_number, group_elem_number,process.shell_string(print_id=False))
373 374 #=========================================================================== 375 # Set the compiler to be gfortran for the loop processes. 376 #===========================================================================
377 - def compiler_choice(self, compiler):
378 """ Different daughter classes might want different compilers. 379 Here, the gfortran compiler is used throughout the compilation 380 (mandatory for CutTools written in f90) """ 381 if isinstance(compiler, str): 382 compiler= {'fortran':compiler, 'f2py':''} 383 384 if not compiler['fortran'] is None and not any([name in compiler['fortran'] for name in \ 385 ['gfortran','ifort']]): 386 logger.info('For loop processes, the compiler must be fortran90'+\ 387 'compatible, like gfortran.') 388 compiler['fortran'] = 'gfortran' 389 self.set_compiler(compiler,True) 390 else: 391 self.set_compiler(compiler)
392
393 - def turn_to_mp_calls(self, helas_calls_list):
394 # Prepend 'MP_' to all the helas calls in helas_calls_list. 395 # Might look like a brutal unsafe implementation, but it is not as 396 # these calls are built from the properties of the HELAS objects and 397 # whether they are evaluated in double or quad precision is none of 398 # their business but only relevant to the output algorithm. 399 # Also the cast to complex masses DCMPLX(*) must be replaced by 400 # CMPLX(*,KIND=16) 401 MP=re.compile(r"(?P<toSub>^.*CALL\s+)",re.IGNORECASE | re.MULTILINE) 402 403 def replaceWith(match_obj): 404 return match_obj.group('toSub')+'MP_'
405 406 DCMPLX=re.compile(r"DCMPLX\((?P<toSub>([^\)]*))\)",\ 407 re.IGNORECASE | re.MULTILINE) 408 409 for i, helas_call in enumerate(helas_calls_list): 410 new_helas_call=MP.sub(replaceWith,helas_call) 411 helas_calls_list[i]=DCMPLX.sub(r"CMPLX(\g<toSub>,KIND=16)",\ 412 new_helas_call)
413 417 425
426 - def make(self):
427 """ Compiles the additional dependences for loop (such as CutTools).""" 428 super(LoopProcessExporterFortranSA, self).make() 429 430 # make CutTools (only necessary with MG option output_dependencies='internal') 431 libdir = os.path.join(self.dir_path,'lib') 432 sourcedir = os.path.join(self.dir_path,'Source') 433 if self.dependencies=='internal': 434 if not os.path.exists(os.path.realpath(pjoin(libdir, 'libcts.a'))) or \ 435 not os.path.exists(os.path.realpath(pjoin(libdir, 'mpmodule.mod'))): 436 if os.path.exists(pjoin(sourcedir,'CutTools')): 437 logger.info('Compiling CutTools (can take a couple of minutes) ...') 438 misc.compile(['CutTools'], cwd = sourcedir) 439 logger.info(' ...done.') 440 else: 441 raise MadGraph5Error('Could not compile CutTools because its'+\ 442 ' source directory could not be found in the SOURCE folder.') 443 if not os.path.exists(os.path.realpath(pjoin(libdir, 'libcts.a'))) or \ 444 not os.path.exists(os.path.realpath(pjoin(libdir, 'mpmodule.mod'))): 445 raise MadGraph5Error('CutTools compilation failed.') 446 447 # Verify compatibility between current compiler and the one which was 448 # used when last compiling CutTools (if specified). 449 compiler_log_path = pjoin(os.path.dirname((os.path.realpath(pjoin( 450 libdir, 'libcts.a')))),'compiler_version.log') 451 if os.path.exists(compiler_log_path): 452 compiler_version_used = open(compiler_log_path,'r').read() 453 if not str(misc.get_gfortran_version(misc.detect_current_compiler(\ 454 pjoin(sourcedir,'make_opts')))) in compiler_version_used: 455 if os.path.exists(pjoin(sourcedir,'CutTools')): 456 logger.info('CutTools was compiled with a different fortran'+\ 457 ' compiler. Re-compiling it now...') 458 misc.compile(['cleanCT'], cwd = sourcedir) 459 misc.compile(['CutTools'], cwd = sourcedir) 460 logger.info(' ...done.') 461 else: 462 raise MadGraph5Error("CutTools installation in %s"\ 463 %os.path.realpath(pjoin(libdir, 'libcts.a'))+\ 464 " seems to have been compiled with a different compiler than"+\ 465 " the one specified in MG5_aMC. Please recompile CutTools.")
466
467 - def cat_coeff(self, ff_number, frac, is_imaginary, Nc_power, Nc_value=3):
468 """Concatenate the coefficient information to reduce it to 469 (fraction, is_imaginary) """ 470 471 total_coeff = ff_number * frac * fractions.Fraction(Nc_value) ** Nc_power 472 473 return (total_coeff, is_imaginary)
474
475 - def get_amp_to_jamp_map(self, col_amps, n_amps):
476 """ Returns a list with element 'i' being a list of tuples corresponding 477 to all apparition of amplitude number 'i' in the jamp number 'j' 478 with coeff 'coeff_j'. The format of each tuple describing an apparition 479 is (j, coeff_j). where coeff_j is of the form (Fraction, is_imag).""" 480 481 if(isinstance(col_amps,list)): 482 if(col_amps and isinstance(col_amps[0],list)): 483 color_amplitudes=col_amps 484 else: 485 raise MadGraph5Error, "Incorrect col_amps argument passed to get_amp_to_jamp_map" 486 else: 487 raise MadGraph5Error, "Incorrect col_amps argument passed to get_amp_to_jamp_map" 488 489 # To store the result 490 res_list = [[] for i in range(n_amps)] 491 for i, coeff_list in enumerate(color_amplitudes): 492 for (coefficient, amp_number) in coeff_list: 493 res_list[amp_number-1].append((i,self.cat_coeff(\ 494 coefficient[0],coefficient[1],coefficient[2],coefficient[3]))) 495 496 return res_list
497
498 - def get_color_matrix(self, matrix_element):
499 """Return the color matrix definition lines. This color matrix is of size 500 NLOOPAMPSxNBORNAMPS and allows for squaring individually each Loop and Born 501 amplitude.""" 502 503 logger.info('Computing diagram color coefficients') 504 505 # The two lists have a list of tuples at element 'i' which correspond 506 # to all apparitions of loop amplitude number 'i' in the jampl number 'j' 507 # with coeff 'coeffj'. The format of each tuple describing an apparition 508 # is (j, coeffj). 509 ampl_to_jampl=self.get_amp_to_jamp_map(\ 510 matrix_element.get_loop_color_amplitudes(), 511 matrix_element.get_number_of_loop_amplitudes()) 512 if matrix_element.get('processes')[0].get('has_born'): 513 ampb_to_jampb=self.get_amp_to_jamp_map(\ 514 matrix_element.get_born_color_amplitudes(), 515 matrix_element.get_number_of_born_amplitudes()) 516 else: 517 ampb_to_jampb=ampl_to_jampl 518 # Below is the original color matrix multiplying the JAMPS 519 if matrix_element.get('color_matrix'): 520 ColorMatrixDenom = \ 521 matrix_element.get('color_matrix').get_line_denominators() 522 ColorMatrixNum = [ matrix_element.get('color_matrix').\ 523 get_line_numerators(index, denominator) for 524 (index, denominator) in enumerate(ColorMatrixDenom) ] 525 else: 526 ColorMatrixDenom= [1] 527 ColorMatrixNum = [[1]] 528 529 # Below is the final color matrix output 530 ColorMatrixNumOutput=[] 531 ColorMatrixDenomOutput=[] 532 533 # Now we construct the color factors between each born and loop amplitude 534 # by scanning their contributions to the different jamps. 535 start = time.time() 536 progress_bar = None 537 time_info = False 538 for i, jampl_list in enumerate(ampl_to_jampl): 539 # This can be pretty long for processes with many color flows. 540 # So, if necessary (i.e. for more than 15s), we tell the user the 541 # estimated time for the processing. 542 if i==5: 543 elapsed_time = time.time()-start 544 t = len(ampl_to_jampl)*(elapsed_time/5.0) 545 if t > 10.0: 546 time_info = True 547 logger.info('The color factors computation will take '+\ 548 ' about %s to run. '%str(datetime.timedelta(seconds=int(t)))+\ 549 'Started on %s.'%datetime.datetime.now().strftime(\ 550 "%d-%m-%Y %H:%M")) 551 if logger.getEffectiveLevel()<logging.WARNING: 552 widgets = ['Color computation:', pbar.Percentage(), ' ', 553 pbar.Bar(),' ', pbar.ETA(), ' '] 554 progress_bar = pbar.ProgressBar(widgets=widgets, 555 maxval=len(ampl_to_jampl), fd=sys.stdout) 556 557 if not progress_bar is None: 558 progress_bar.update(i+1) 559 # Flush to force the printout of the progress_bar to be updated 560 sys.stdout.flush() 561 562 line_num=[] 563 line_denom=[] 564 565 # Treat the special case where this specific amplitude contributes to no 566 # color flow at all. So it is zero because of color but not even due to 567 # an accidental cancellation among color flows, but simply because of its 568 # projection to each individual color flow is zero. In such case, the 569 # corresponding jampl_list is empty and all color coefficients must then 570 # be zero. This happens for example in the Higgs Effective Theory model 571 # for the bubble made of a 4-gluon vertex and the effective ggH vertex. 572 if len(jampl_list)==0: 573 line_num=[0]*len(ampb_to_jampb) 574 line_denom=[1]*len(ampb_to_jampb) 575 ColorMatrixNumOutput.append(line_num) 576 ColorMatrixDenomOutput.append(line_denom) 577 continue 578 579 for jampb_list in ampb_to_jampb: 580 real_num=0 581 imag_num=0 582 common_denom=color_amp.ColorMatrix.lcmm(*[abs(ColorMatrixDenom[jampl]* 583 ampl_coeff[0].denominator*ampb_coeff[0].denominator) for 584 ((jampl, ampl_coeff),(jampb,ampb_coeff)) in 585 itertools.product(jampl_list,jampb_list)]) 586 for ((jampl, ampl_coeff),(jampb, ampb_coeff)) in \ 587 itertools.product(jampl_list,jampb_list): 588 # take the numerator and multiply by lcm/denominator 589 # as we will later divide by the lcm. 590 buff_num=ampl_coeff[0].numerator*\ 591 ampb_coeff[0].numerator*ColorMatrixNum[jampl][jampb]*\ 592 abs(common_denom)/(ampl_coeff[0].denominator*\ 593 ampb_coeff[0].denominator*ColorMatrixDenom[jampl]) 594 # Remember that we must take the complex conjugate of 595 # the born jamp color coefficient because we will compute 596 # the square with 2 Re(LoopAmp x BornAmp*) 597 if ampl_coeff[1] and ampb_coeff[1]: 598 real_num=real_num+buff_num 599 elif not ampl_coeff[1] and not ampb_coeff[1]: 600 real_num=real_num+buff_num 601 elif not ampl_coeff[1] and ampb_coeff[1]: 602 imag_num=imag_num-buff_num 603 else: 604 imag_num=imag_num+buff_num 605 assert not (real_num!=0 and imag_num!=0), "MadGraph5_aMC@NLO found a "+\ 606 "color matrix element which has both a real and imaginary part." 607 if imag_num!=0: 608 res=fractions.Fraction(imag_num,common_denom) 609 line_num.append(res.numerator) 610 # Negative denominator means imaginary color coef of the 611 # final color matrix 612 line_denom.append(res.denominator*-1) 613 else: 614 res=fractions.Fraction(real_num,common_denom) 615 line_num.append(res.numerator) 616 # Positive denominator means real color coef of the final color matrix 617 line_denom.append(res.denominator) 618 619 ColorMatrixNumOutput.append(line_num) 620 ColorMatrixDenomOutput.append(line_denom) 621 622 if time_info: 623 logger.info('Finished on %s.'%datetime.datetime.now().strftime(\ 624 "%d-%m-%Y %H:%M")) 625 if progress_bar!=None: 626 progress_bar.finish() 627 628 return (ColorMatrixNumOutput,ColorMatrixDenomOutput)
629
630 - def get_context(self,matrix_element):
631 """ Returns the contextual variables which need to be set when 632 pre-processing the template files.""" 633 634 # The nSquaredSO entry of the general replace dictionary should have 635 # been set in write_loopmatrix prior to the first call to this function 636 # However, for cases where the TIRCaching contextual variable is 637 # irrelevant (like in the default output), this might not be the case 638 # so we set it to 1. 639 try: 640 n_squared_split_orders = matrix_element.rep_dict['nSquaredSO'] 641 except (KeyError, AttributeError): 642 n_squared_split_orders = 1 643 644 LoopInduced = not matrix_element.get('processes')[0].get('has_born') 645 646 # Force the computation of loop color flows for loop_induced processes 647 ComputeColorFlows = self.compute_color_flows or LoopInduced 648 # The variable AmplitudeReduction is just to make the contextual 649 # conditions more readable in the include files. 650 AmplitudeReduction = LoopInduced or ComputeColorFlows 651 # Even when not reducing at the amplitude level, the TIR caching 652 # is useful when there is more than one squared split order config. 653 TIRCaching = AmplitudeReduction or n_squared_split_orders>1 654 MadEventOutput = False 655 656 return {'LoopInduced': LoopInduced, 657 'ComputeColorFlows': ComputeColorFlows, 658 'AmplitudeReduction': AmplitudeReduction, 659 'TIRCaching': TIRCaching, 660 'MadEventOutput': MadEventOutput}
661 662 #=========================================================================== 663 # generate_subprocess_directory_v4 664 #===========================================================================
665 - def generate_loop_subprocess(self, matrix_element, fortran_model, 666 group_number = None, proc_id = None, config_map=None):
667 """Generate the Pxxxxx directory for a loop subprocess in MG4 standalone, 668 including the necessary loop_matrix.f, born_matrix.f and include files. 669 Notice that this is too different from generate_subprocess_directory_v4 670 so that there is no point reusing this mother function. 671 The 'group_number' and 'proc_id' options are only used for the LoopInduced 672 MadEvent output and only to specify the ME_identifier and the P* 673 SubProcess directory name.""" 674 675 cwd = os.getcwd() 676 proc_dir_name = self.get_SubProc_folder_name( 677 matrix_element.get('processes')[0],group_number,proc_id) 678 dirpath = os.path.join(self.dir_path, 'SubProcesses', proc_dir_name) 679 680 try: 681 os.mkdir(dirpath) 682 except os.error as error: 683 logger.warning(error.strerror + " " + dirpath) 684 685 try: 686 os.chdir(dirpath) 687 except os.error: 688 logger.error('Could not cd to directory %s' % dirpath) 689 return 0 690 691 logger.info('Creating files in directory %s' % dirpath) 692 693 # Extract number of external particles 694 (nexternal, ninitial) = matrix_element.get_nexternal_ninitial() 695 696 calls=self.write_loop_matrix_element_v4(None,matrix_element, 697 fortran_model, group_number = group_number, 698 proc_id = proc_id, config_map = config_map) 699 700 # We assume here that all processes must share the same property of 701 # having a born or not, which must be true anyway since these are two 702 # definite different classes of processes which can never be treated on 703 # the same footing. 704 if matrix_element.get('processes')[0].get('has_born'): 705 filename = 'born_matrix.f' 706 calls = self.write_bornmatrix( 707 writers.FortranWriter(filename), 708 matrix_element, 709 fortran_model) 710 711 filename = 'pmass.inc' 712 self.write_pmass_file(writers.FortranWriter(filename), 713 matrix_element) 714 715 filename = 'ngraphs.inc' 716 self.write_ngraphs_file(writers.FortranWriter(filename), 717 len(matrix_element.get_all_amplitudes())) 718 719 # Do not draw the loop diagrams if they are too many. 720 # The user can always decide to do it manually, if really needed 721 loop_diags = [loop_diag for loop_diag in\ 722 matrix_element.get('base_amplitude').get('loop_diagrams')\ 723 if isinstance(loop_diag,LoopDiagram) and loop_diag.get('type') > 0] 724 if len(loop_diags)>5000: 725 logger.info("There are more than 5000 loop diagrams."+\ 726 "Only the first 5000 are drawn.") 727 filename = "loop_matrix.ps" 728 plot = draw.MultiEpsDiagramDrawer(base_objects.DiagramList( 729 loop_diags[:5000]),filename, 730 model=matrix_element.get('processes')[0].get('model'),amplitude='') 731 logger.info("Drawing loop Feynman diagrams for " + \ 732 matrix_element.get('processes')[0].nice_string()) 733 plot.draw() 734 735 if matrix_element.get('processes')[0].get('has_born'): 736 filename = "born_matrix.ps" 737 plot = draw.MultiEpsDiagramDrawer(matrix_element.get('base_amplitude').\ 738 get('born_diagrams'), 739 filename, 740 model=matrix_element.get('processes')[0].\ 741 get('model'), 742 amplitude='') 743 logger.info("Generating born Feynman diagrams for " + \ 744 matrix_element.get('processes')[0].nice_string(\ 745 print_weighted=False)) 746 plot.draw() 747 748 self.link_files_from_Subprocesses(self.get_SubProc_folder_name( 749 matrix_element.get('processes')[0],group_number,proc_id)) 750 751 # Return to original PWD 752 os.chdir(cwd) 753 754 if not calls: 755 calls = 0 756 return calls
757 778
779 - def generate_general_replace_dict(self,matrix_element, 780 group_number = None, proc_id = None):
781 """Generates the entries for the general replacement dictionary used 782 for the different output codes for this exporter.The arguments 783 group_number and proc_id are just for the LoopInduced output with MadEvent.""" 784 785 dict={} 786 # A general process prefix which appears in front of all MadLooop 787 # subroutines and common block so that several processes can be compiled 788 # together into one library, as necessary to follow BLHA guidelines. 789 790 dict['proc_prefix'] = self.get_ME_identifier(matrix_element, 791 group_number = group_number, group_elem_number = proc_id) 792 793 # The proc_id is used for MadEvent grouping, so none of our concern here 794 # and it is simply set to an empty string. 795 dict['proc_id'] = '' 796 # Extract version number and date from VERSION file 797 info_lines = self.get_mg5_info_lines() 798 dict['info_lines'] = info_lines 799 # Extract process info lines 800 process_lines = self.get_process_info_lines(matrix_element) 801 dict['process_lines'] = process_lines 802 # Extract number of external particles 803 (nexternal, ninitial) = matrix_element.get_nexternal_ninitial() 804 dict['nexternal'] = nexternal 805 dict['nincoming'] = ninitial 806 # Extract ncomb 807 ncomb = matrix_element.get_helicity_combinations() 808 dict['ncomb'] = ncomb 809 # Extract nloopamps 810 nloopamps = matrix_element.get_number_of_loop_amplitudes() 811 dict['nloopamps'] = nloopamps 812 # Extract nloopdiags 813 nloopdiags = len(matrix_element.get('diagrams')) 814 dict['nloopdiags'] = nloopdiags 815 # Extract nctamps 816 nctamps = matrix_element.get_number_of_CT_amplitudes() 817 dict['nctamps'] = nctamps 818 # Extract nwavefuncs 819 nwavefuncs = matrix_element.get_number_of_external_wavefunctions() 820 dict['nwavefuncs'] = nwavefuncs 821 # Set format of the double precision 822 dict['real_dp_format']='real*8' 823 dict['real_mp_format']='real*16' 824 # Set format of the complex 825 dict['complex_dp_format']='complex*16' 826 dict['complex_mp_format']='complex*32' 827 # Set format of the masses 828 dict['mass_dp_format'] = dict['complex_dp_format'] 829 dict['mass_mp_format'] = dict['complex_mp_format'] 830 # Fill in default values for the placeholders for the madevent 831 # loop-induced output 832 dict['nmultichannels'] = 0 833 dict['nmultichannel_configs'] = 0 834 dict['config_map_definition'] = '' 835 dict['config_index_map_definition'] = '' 836 # Color matrix size 837 # For loop induced processes it is NLOOPAMPSxNLOOPAMPS and otherwise 838 # it is NLOOPAMPSxNBORNAMPS 839 # Also, how to access the number of Born squared order contributions 840 841 if matrix_element.get('processes')[0].get('has_born'): 842 dict['color_matrix_size'] = 'nbornamps' 843 dict['get_nsqso_born']=\ 844 "include 'nsqso_born.inc'" 845 else: 846 dict['get_nsqso_born']="""INTEGER NSQSO_BORN 847 PARAMETER (NSQSO_BORN=0) 848 """ 849 dict['color_matrix_size'] = 'nloopamps' 850 851 # These placeholders help to have as many common templates for the 852 # output of the loop induced processes and those with a born 853 # contribution. 854 if matrix_element.get('processes')[0].get('has_born'): 855 # Extract nbornamps 856 nbornamps = matrix_element.get_number_of_born_amplitudes() 857 dict['nbornamps'] = nbornamps 858 dict['ncomb_helas_objs'] = ',ncomb' 859 dict['nbornamps_decl'] = \ 860 """INTEGER NBORNAMPS 861 PARAMETER (NBORNAMPS=%d)"""%nbornamps 862 dict['nBornAmps'] = nbornamps 863 864 else: 865 dict['ncomb_helas_objs'] = '' 866 dict['dp_born_amps_decl'] = '' 867 dict['dp_born_amps_decl_in_mp'] = '' 868 dict['copy_mp_to_dp_born_amps'] = '' 869 dict['mp_born_amps_decl'] = '' 870 dict['nbornamps_decl'] = '' 871 dict['nbornamps'] = 0 872 dict['nBornAmps'] = 0 873 874 return dict
875
876 - def write_loop_matrix_element_v4(self, writer, matrix_element, fortran_model, 877 group_number = None, proc_id = None, config_map = None):
878 """ Writes loop_matrix.f, CT_interface.f, loop_num.f and 879 mp_born_amps_and_wfs. 880 The arguments group_number and proc_id are just for the LoopInduced 881 output with MadEvent and only used in get_ME_identifier.""" 882 883 # Create the necessary files for the loop matrix element subroutine 884 885 if config_map: 886 raise MadGraph5Error, 'The default loop output cannot be used with'+\ 887 'MadEvent and cannot compute the AMP2 for multi-channeling.' 888 889 if not isinstance(fortran_model,\ 890 helas_call_writers.FortranUFOHelasCallWriter): 891 raise MadGraph5Error, 'The loop fortran output can only'+\ 892 ' work with a UFO Fortran model' 893 894 LoopFortranModel = helas_call_writers.FortranUFOHelasCallWriter( 895 argument=fortran_model.get('model'), 896 hel_sum=matrix_element.get('processes')[0].get('has_born')) 897 898 # Compute the analytical information of the loop wavefunctions in the 899 # loop helas matrix elements using the cached aloha model to reuse 900 # as much as possible the aloha computations already performed for 901 # writing out the aloha fortran subroutines. 902 matrix_element.compute_all_analytic_information( 903 self.get_aloha_model(matrix_element.get('processes')[0].get('model'))) 904 905 # Initialize a general replacement dictionary with entries common to 906 # many files generated here. 907 matrix_element.rep_dict = self.generate_general_replace_dict( 908 matrix_element, group_number = group_number, proc_id = proc_id) 909 910 # Extract max number of loop couplings (specific to this output type) 911 matrix_element.rep_dict['maxlcouplings']= \ 912 matrix_element.find_max_loop_coupling() 913 # The born amp declaration suited for also outputing the loop-induced 914 # processes as well. 915 if matrix_element.get('processes')[0].get('has_born'): 916 matrix_element.rep_dict['dp_born_amps_decl_in_mp'] = \ 917 matrix_element.rep_dict['complex_dp_format']+" DPAMP(NBORNAMPS,NCOMB)"+\ 918 "\n common/%sAMPS/DPAMP"%matrix_element.rep_dict['proc_prefix'] 919 matrix_element.rep_dict['dp_born_amps_decl'] = \ 920 matrix_element.rep_dict['complex_dp_format']+" AMP(NBORNAMPS,NCOMB)"+\ 921 "\n common/%sAMPS/AMP"%matrix_element.rep_dict['proc_prefix'] 922 matrix_element.rep_dict['mp_born_amps_decl'] = \ 923 matrix_element.rep_dict['complex_mp_format']+" AMP(NBORNAMPS,NCOMB)"+\ 924 "\n common/%sMP_AMPS/AMP"%matrix_element.rep_dict['proc_prefix'] 925 matrix_element.rep_dict['copy_mp_to_dp_born_amps'] = \ 926 '\n'.join(['DO I=1,NBORNAMPS','DPAMP(I,H)=AMP(I,H)','ENDDO']) 927 928 if writer: 929 raise MadGraph5Error, 'Matrix output mode no longer supported.' 930 931 filename = 'loop_matrix.f' 932 calls = self.write_loopmatrix(writers.FortranWriter(filename), 933 matrix_element, 934 LoopFortranModel) 935 936 # Write out the proc_prefix in a file, this is quite handy 937 proc_prefix_writer = writers.FortranWriter('proc_prefix.txt','w') 938 proc_prefix_writer.write(matrix_element.rep_dict['proc_prefix']) 939 proc_prefix_writer.close() 940 941 filename = 'check_sa.f' 942 self.write_check_sa(writers.FortranWriter(filename),matrix_element) 943 944 filename = 'CT_interface.f' 945 self.write_CT_interface(writers.FortranWriter(filename),\ 946 matrix_element) 947 948 filename = 'improve_ps.f' 949 calls = self.write_improve_ps(writers.FortranWriter(filename), 950 matrix_element) 951 952 filename = 'loop_num.f' 953 self.write_loop_num(writers.FortranWriter(filename),\ 954 matrix_element,LoopFortranModel) 955 956 filename = 'mp_born_amps_and_wfs.f' 957 self.write_born_amps_and_wfs(writers.FortranWriter(filename),\ 958 matrix_element,LoopFortranModel) 959 960 # Extract number of external particles 961 (nexternal, ninitial) = matrix_element.get_nexternal_ninitial() 962 filename = 'nexternal.inc' 963 self.write_nexternal_file(writers.FortranWriter(filename), 964 nexternal, ninitial) 965 966 return calls
967
968 - def generate_subprocess_directory_v4(self, matrix_element, fortran_model):
969 """ To overload the default name for this function such that the correct 970 function is used when called from the command interface """ 971 972 return self.generate_loop_subprocess(matrix_element,fortran_model)
973
974 - def write_check_sa(self, writer, matrix_element):
975 """Writes out the steering code check_sa. In the optimized output mode, 976 All the necessary entries in the replace_dictionary have already been 977 set in write_loopmatrix because it is only there that one has access to 978 the information about split orders.""" 979 replace_dict = copy.copy(matrix_element.rep_dict) 980 for key in ['print_so_born_results','print_so_loop_results', 981 'write_so_born_results','write_so_loop_results','set_coupling_target']: 982 if key not in replace_dict.keys(): 983 replace_dict[key]='' 984 if matrix_element.get('processes')[0].get('has_born'): 985 file = open(os.path.join(self.template_dir,'check_sa.inc')).read() 986 elif self.opt['mode'] == 'reweight': 987 file = open(os.path.join(self.template_dir,\ 988 'check_py.f')).read() 989 else: 990 file = open(os.path.join(self.template_dir,\ 991 'check_sa_loop_induced.inc')).read() 992 file=file%replace_dict 993 writer.writelines(file)
994
995 - def write_improve_ps(self, writer, matrix_element):
996 """ Write out the improve_ps subroutines which modify the PS point 997 given in input and slightly deform it to achieve exact onshellness on 998 all external particles as well as perfect energy-momentum conservation""" 999 replace_dict = copy.copy(matrix_element.rep_dict) 1000 1001 (nexternal,ninitial)=matrix_element.get_nexternal_ninitial() 1002 replace_dict['ninitial']=ninitial 1003 mass_list=matrix_element.get_external_masses()[:-2] 1004 mp_variable_prefix = check_param_card.ParamCard.mp_prefix 1005 1006 # Write the quadruple precision version of this routine only. 1007 replace_dict['real_format']=replace_dict['real_mp_format'] 1008 replace_dict['mp_prefix']='MP_' 1009 replace_dict['exp_letter']='e' 1010 replace_dict['mp_specifier']='_16' 1011 replace_dict['coupl_inc_name']='mp_coupl.inc' 1012 replace_dict['masses_def']='\n'.join(['MASSES(%(i)d)=%(prefix)s%(m)s'\ 1013 %{'i':i+1,'m':m, 'prefix':mp_variable_prefix} for \ 1014 i, m in enumerate(mass_list)]) 1015 file_mp = open(os.path.join(self.template_dir,'improve_ps.inc')).read() 1016 file_mp=file_mp%replace_dict 1017 # 1018 writer.writelines(file_mp)
1019
1020 - def write_loop_num(self, writer, matrix_element,fortran_model):
1021 """ Create the file containing the core subroutine called by CutTools 1022 which contains the Helas calls building the loop""" 1023 1024 if not matrix_element.get('processes') or \ 1025 not matrix_element.get('diagrams'): 1026 return 0 1027 1028 # Set lowercase/uppercase Fortran code 1029 writers.FortranWriter.downcase = False 1030 1031 file = open(os.path.join(self.template_dir,'loop_num.inc')).read() 1032 1033 replace_dict = copy.copy(matrix_element.rep_dict) 1034 1035 loop_helas_calls=fortran_model.get_loop_amplitude_helas_calls(matrix_element) 1036 replace_dict['maxlcouplings']=matrix_element.find_max_loop_coupling() 1037 replace_dict['loop_helas_calls'] = "\n".join(loop_helas_calls) 1038 1039 # The squaring is only necessary for the processes with born where the 1040 # sum over helicities is done before sending the numerator to CT. 1041 dp_squaring_lines=['DO I=1,NBORNAMPS', 1042 'CFTOT=DCMPLX(CF_N(AMPLNUM,I)/DBLE(ABS(CF_D(AMPLNUM,I))),0.0d0)', 1043 'IF(CF_D(AMPLNUM,I).LT.0) CFTOT=CFTOT*IMAG1', 1044 'RES=RES+CFTOT*BUFF*DCONJG(AMP(I,H))','ENDDO'] 1045 mp_squaring_lines=['DO I=1,NBORNAMPS', 1046 'CFTOT=CMPLX(CF_N(AMPLNUM,I)/(1.0E0_16*ABS(CF_D(AMPLNUM,I))),0.0E0_16,KIND=16)', 1047 'IF(CF_D(AMPLNUM,I).LT.0) CFTOT=CFTOT*IMAG1', 1048 'QPRES=QPRES+CFTOT*BUFF*CONJG(AMP(I,H))','ENDDO'] 1049 if matrix_element.get('processes')[0].get('has_born'): 1050 replace_dict['dp_squaring']='\n'.join(dp_squaring_lines) 1051 replace_dict['mp_squaring']='\n'.join(mp_squaring_lines) 1052 else: 1053 replace_dict['dp_squaring']='RES=BUFF' 1054 replace_dict['mp_squaring']='QPRES=BUFF' 1055 1056 # Prepend MP_ to all helas calls. 1057 self.turn_to_mp_calls(loop_helas_calls) 1058 replace_dict['mp_loop_helas_calls'] = "\n".join(loop_helas_calls) 1059 1060 file=file%replace_dict 1061 1062 if writer: 1063 writer.writelines(file) 1064 else: 1065 return file
1066
1067 - def write_CT_interface(self, writer, matrix_element, optimized_output=False):
1068 """ Create the file CT_interface.f which contains the subroutine defining 1069 the loop HELAS-like calls along with the general interfacing subroutine. """ 1070 1071 files=[] 1072 1073 # First write CT_interface which interfaces MG5 with CutTools. 1074 replace_dict=copy.copy(matrix_element.rep_dict) 1075 1076 # We finalize CT result differently wether we used the built-in 1077 # squaring against the born. 1078 if matrix_element.get('processes')[0].get('has_born'): 1079 replace_dict['finalize_CT']='\n'.join([\ 1080 'RES(%d)=NORMALIZATION*2.0d0*DBLE(RES(%d))'%(i,i) for i in range(1,4)]) 1081 else: 1082 replace_dict['finalize_CT']='\n'.join([\ 1083 'RES(%d)=NORMALIZATION*RES(%d)'%(i,i) for i in range(1,4)]) 1084 1085 file = open(os.path.join(self.template_dir,'CT_interface.inc')).read() 1086 1087 file = file % replace_dict 1088 files.append(file) 1089 1090 # Now collect the different kind of subroutines needed for the 1091 # loop HELAS-like calls. 1092 HelasLoopAmpsCallKeys=matrix_element.get_used_helas_loop_amps() 1093 1094 for callkey in HelasLoopAmpsCallKeys: 1095 replace_dict=copy.copy(matrix_element.rep_dict) 1096 # Add to this dictionary all other attribute common to all 1097 # HELAS-like loop subroutines. 1098 if matrix_element.get('processes')[0].get('has_born'): 1099 replace_dict['validh_or_nothing']=',validh' 1100 else: 1101 replace_dict['validh_or_nothing']='' 1102 # In the optimized output, the number of couplings in the loop is 1103 # not specified so we only treat it here if necessary: 1104 if len(callkey)>2: 1105 replace_dict['ncplsargs']=callkey[2] 1106 cplsargs="".join(["C%d,MP_C%d, "%(i,i) for i in range(1,callkey[2]+1)]) 1107 replace_dict['cplsargs']=cplsargs 1108 cplsdecl="".join(["C%d, "%i for i in range(1,callkey[2]+1)])[:-2] 1109 replace_dict['cplsdecl']=cplsdecl 1110 mp_cplsdecl="".join(["MP_C%d, "%i for i in range(1,callkey[2]+1)])[:-2] 1111 replace_dict['mp_cplsdecl']=mp_cplsdecl 1112 cplset="\n".join(["\n".join(["LC(%d)=C%d"%(i,i),\ 1113 "MP_LC(%d)=MP_C%d"%(i,i)])\ 1114 for i in range(1,callkey[2]+1)]) 1115 replace_dict['cplset']=cplset 1116 1117 replace_dict['nloopline']=callkey[0] 1118 wfsargs="".join(["W%d, "%i for i in range(1,callkey[1]+1)]) 1119 replace_dict['wfsargs']=wfsargs 1120 # We don't pass the multiple precision mass in the optimized_output 1121 if not optimized_output: 1122 margs="".join(["M%d,MP_M%d, "%(i,i) for i in range(1,callkey[0]+1)]) 1123 else: 1124 margs="".join(["M%d, "%i for i in range(1,callkey[0]+1)]) 1125 replace_dict['margs']=margs 1126 wfsargsdecl="".join([("W%d, "%i) for i in range(1,callkey[1]+1)])[:-2] 1127 replace_dict['wfsargsdecl']=wfsargsdecl 1128 margsdecl="".join(["M%d, "%i for i in range(1,callkey[0]+1)])[:-2] 1129 replace_dict['margsdecl']=margsdecl 1130 mp_margsdecl="".join(["MP_M%d, "%i for i in range(1,callkey[0]+1)])[:-2] 1131 replace_dict['mp_margsdecl']=mp_margsdecl 1132 weset="\n".join([("WE("+str(i)+")=W"+str(i)) for \ 1133 i in range(1,callkey[1]+1)]) 1134 replace_dict['weset']=weset 1135 weset="\n".join([("WE(%d)=W%d"%(i,i)) for i in range(1,callkey[1]+1)]) 1136 replace_dict['weset']=weset 1137 msetlines=["M2L(1)=M%d**2"%(callkey[0]),] 1138 mset="\n".join(msetlines+["M2L(%d)=M%d**2"%(i,i-1) for \ 1139 i in range(2,callkey[0]+1)]) 1140 replace_dict['mset']=mset 1141 mset2lines=["ML(1)=M%d"%(callkey[0]),"ML(2)=M%d"%(callkey[0]), 1142 "MP_ML(1)=MP_M%d"%(callkey[0]),"MP_ML(2)=MP_M%d"%(callkey[0])] 1143 mset2="\n".join(mset2lines+["\n".join(["ML(%d)=M%d"%(i,i-2), 1144 "MP_ML(%d)=MP_M%d"%(i,i-2)]) for \ 1145 i in range(3,callkey[0]+3)]) 1146 replace_dict['mset2']=mset2 1147 replace_dict['nwfsargs'] = callkey[1] 1148 if callkey[0]==callkey[1]: 1149 replace_dict['nwfsargs_header'] = "" 1150 replace_dict['pairingargs']="" 1151 replace_dict['pairingdecl']="" 1152 pairingset="""DO I=1,NLOOPLINE 1153 PAIRING(I)=1 1154 ENDDO 1155 """ 1156 replace_dict['pairingset']=pairingset 1157 else: 1158 replace_dict['nwfsargs_header'] = '_%d'%callkey[1] 1159 pairingargs="".join([("P"+str(i)+", ") for i in \ 1160 range(1,callkey[0]+1)]) 1161 replace_dict['pairingargs']=pairingargs 1162 pairingdecl="integer "+"".join([("P"+str(i)+", ") for i in \ 1163 range(1,callkey[0]+1)])[:-2] 1164 replace_dict['pairingdecl']=pairingdecl 1165 pairingset="\n".join([("PAIRING("+str(i)+")=P"+str(i)) for \ 1166 i in range(1,callkey[0]+1)]) 1167 replace_dict['pairingset']=pairingset 1168 1169 file = open(os.path.join(self.template_dir,\ 1170 'helas_loop_amplitude.inc')).read() 1171 file = file % replace_dict 1172 files.append(file) 1173 1174 file="\n".join(files) 1175 1176 if writer: 1177 writer.writelines(file,context=self.get_context(matrix_element)) 1178 else: 1179 return file
1180 1181 # Helper function to split HELAS CALLS in dedicated subroutines placed 1182 # in different files.
1183 - def split_HELASCALLS(self, writer, replace_dict, template_name, masterfile, \ 1184 helas_calls, entry_name, bunch_name,n_helas=2000, 1185 required_so_broadcaster = 'LOOP_REQ_SO_DONE', 1186 continue_label = 1000, context={}):
1187 """ Finish the code generation with splitting. 1188 Split the helas calls in the argument helas_calls into bunches of 1189 size n_helas and place them in dedicated subroutine with name 1190 <bunch_name>_i. Also setup the corresponding calls to these subroutine 1191 in the replace_dict dictionary under the entry entry_name. 1192 The context specified will be forwarded to the the fileWriter.""" 1193 helascalls_replace_dict=copy.copy(replace_dict) 1194 helascalls_replace_dict['bunch_name']=bunch_name 1195 helascalls_files=[] 1196 for i, k in enumerate(range(0, len(helas_calls), n_helas)): 1197 helascalls_replace_dict['bunch_number']=i+1 1198 helascalls_replace_dict['helas_calls']=\ 1199 '\n'.join(helas_calls[k:k + n_helas]) 1200 helascalls_replace_dict['required_so_broadcaster']=\ 1201 required_so_broadcaster 1202 helascalls_replace_dict['continue_label']=continue_label 1203 new_helascalls_file = open(os.path.join(self.template_dir,\ 1204 template_name)).read() 1205 new_helascalls_file = new_helascalls_file % helascalls_replace_dict 1206 helascalls_files.append(new_helascalls_file) 1207 # Setup the call to these HELASCALLS subroutines in loop_matrix.f 1208 helascalls_calls = [ "CALL %s%s_%d(P,NHEL,H,IC)"%\ 1209 (replace_dict['proc_prefix'] ,bunch_name,a+1) \ 1210 for a in range(len(helascalls_files))] 1211 replace_dict[entry_name]='\n'.join(helascalls_calls) 1212 if writer: 1213 for i, helascalls_file in enumerate(helascalls_files): 1214 filename = '%s_%d.f'%(bunch_name,i+1) 1215 writers.FortranWriter(filename).writelines(helascalls_file, 1216 context=context) 1217 else: 1218 masterfile='\n'.join([masterfile,]+helascalls_files) 1219 1220 return masterfile
1221
1222 - def write_loopmatrix(self, writer, matrix_element, fortran_model, \ 1223 noSplit=False):
1224 """Create the loop_matrix.f file.""" 1225 1226 if not matrix_element.get('processes') or \ 1227 not matrix_element.get('diagrams'): 1228 return 0 1229 1230 # Set lowercase/uppercase Fortran code 1231 1232 writers.FortranWriter.downcase = False 1233 1234 replace_dict = copy.copy(matrix_element.rep_dict) 1235 1236 # Extract overall denominator 1237 # Averaging initial state color, spin, and identical FS particles 1238 den_factor_line = self.get_den_factor_line(matrix_element) 1239 replace_dict['den_factor_line'] = den_factor_line 1240 # When the user asks for the polarized matrix element we must 1241 # multiply back by the helicity averaging factor 1242 replace_dict['hel_avg_factor'] = matrix_element.get_hel_avg_factor() 1243 1244 # These entries are specific for the output for loop-induced processes 1245 # Also sets here the details of the squaring of the loop ampltiudes 1246 # with the born or the loop ones. 1247 if not matrix_element.get('processes')[0].get('has_born'): 1248 replace_dict['compute_born']=\ 1249 """C There is of course no born for loop induced processes 1250 ANS(0)=0.0d0 1251 """ 1252 replace_dict['set_reference']='\n'.join([ 1253 'C For loop-induced, the reference for comparison is set later'+\ 1254 ' from the total contribution of the previous PS point considered.', 1255 'C But you can edit here the value to be used for the first PS point.', 1256 'if (NPSPOINTS.eq.0) then','ref=1.0d-50','else', 1257 'ref=nextRef/DBLE(NPSPOINTS)','endif']) 1258 replace_dict['loop_induced_setup'] = '\n'.join([ 1259 'HELPICKED_BU=HELPICKED','HELPICKED=H','MP_DONE=.FALSE.', 1260 'IF(SKIPLOOPEVAL) THEN','GOTO 1227','ENDIF']) 1261 replace_dict['loop_induced_finalize'] = \ 1262 ("""DO I=NCTAMPS+1,NLOOPAMPS 1263 IF((CTMODERUN.NE.-1).AND..NOT.CHECKPHASE.AND.(.NOT.S(I))) THEN 1264 WRITE(*,*) '##W03 WARNING Contribution ',I 1265 WRITE(*,*) ' is unstable for helicity ',H 1266 ENDIF 1267 C IF(.NOT.%(proc_prefix)sISZERO(ABS(AMPL(2,I))+ABS(AMPL(3,I)),REF,-1,H)) THEN 1268 C WRITE(*,*) '##W04 WARNING Contribution ',I,' for helicity ',H,' has a contribution to the poles.' 1269 C WRITE(*,*) 'Finite contribution = ',AMPL(1,I) 1270 C WRITE(*,*) 'single pole contribution = ',AMPL(2,I) 1271 C WRITE(*,*) 'double pole contribution = ',AMPL(3,I) 1272 C ENDIF 1273 ENDDO 1274 1227 CONTINUE 1275 HELPICKED=HELPICKED_BU""")%replace_dict 1276 replace_dict['loop_helas_calls']="" 1277 replace_dict['nctamps_or_nloopamps']='nloopamps' 1278 replace_dict['nbornamps_or_nloopamps']='nloopamps' 1279 replace_dict['squaring']=\ 1280 """ANS(1)=ANS(1)+DBLE(CFTOT*AMPL(1,I)*DCONJG(AMPL(1,J))) 1281 IF (J.EQ.1) THEN 1282 ANS(2)=ANS(2)+DBLE(CFTOT*AMPL(2,I))+DIMAG(CFTOT*AMPL(2,I)) 1283 ANS(3)=ANS(3)+DBLE(CFTOT*AMPL(3,I))+DIMAG(CFTOT*AMPL(3,I)) 1284 ENDIF""" 1285 else: 1286 replace_dict['compute_born']=\ 1287 """C Compute the born, for a specific helicity if asked so. 1288 call %(proc_prefix)ssmatrixhel(P_USER,USERHEL,ANS(0)) 1289 """%matrix_element.rep_dict 1290 replace_dict['set_reference']=\ 1291 """C We chose to use the born evaluation for the reference 1292 call %(proc_prefix)ssmatrix(p,ref)"""%matrix_element.rep_dict 1293 replace_dict['loop_induced_helas_calls'] = "" 1294 replace_dict['loop_induced_finalize'] = "" 1295 replace_dict['loop_induced_setup'] = "" 1296 replace_dict['nctamps_or_nloopamps']='nctamps' 1297 replace_dict['nbornamps_or_nloopamps']='nbornamps' 1298 replace_dict['squaring']='\n'.join(['DO K=1,3', 1299 'ANS(K)=ANS(K)+2.0d0*DBLE(CFTOT*AMPL(K,I)*DCONJG(AMP(J,H)))', 1300 'ENDDO']) 1301 1302 # Write a dummy nsquaredSO.inc which is used in the default 1303 # loop_matrix.f code (even though it does not support split orders evals) 1304 # just to comply with the syntax expected from the external code using MadLoop. 1305 writers.FortranWriter('nsquaredSO.inc').writelines( 1306 """INTEGER NSQUAREDSO 1307 PARAMETER (NSQUAREDSO=0)""") 1308 1309 # Actualize results from the loops computed. Only necessary for 1310 # processes with a born. 1311 actualize_ans=[] 1312 if matrix_element.get('processes')[0].get('has_born'): 1313 actualize_ans.append("DO I=NCTAMPS+1,NLOOPAMPS") 1314 actualize_ans.extend("ANS(%d)=ANS(%d)+AMPL(%d,I)"%(i,i,i) for i \ 1315 in range(1,4)) 1316 actualize_ans.append(\ 1317 "IF((CTMODERUN.NE.-1).AND..NOT.CHECKPHASE.AND.(.NOT.S(I))) THEN") 1318 actualize_ans.append(\ 1319 "WRITE(*,*) '##W03 WARNING Contribution ',I,' is unstable.'") 1320 actualize_ans.extend(["ENDIF","ENDDO"]) 1321 replace_dict['actualize_ans']='\n'.join(actualize_ans) 1322 else: 1323 replace_dict['actualize_ans']=\ 1324 ("""C We add five powers to the reference value to loosen a bit the vanishing pole check. 1325 C IF(.NOT.(CHECKPHASE.OR.(.NOT.HELDOUBLECHECKED)).AND..NOT.%(proc_prefix)sISZERO(ABS(ANS(2))+ABS(ANS(3)),ABS(ANS(1))*(10.0d0**5),-1,H)) THEN 1326 C WRITE(*,*) '##W05 WARNING Found a PS point with a contribution to the single pole.' 1327 C WRITE(*,*) 'Finite contribution = ',ANS(1) 1328 C WRITE(*,*) 'single pole contribution = ',ANS(2) 1329 C WRITE(*,*) 'double pole contribution = ',ANS(3) 1330 C ENDIF""")%replace_dict 1331 1332 # Write out the color matrix 1333 (CMNum,CMDenom) = self.get_color_matrix(matrix_element) 1334 CMWriter=open(pjoin('..','MadLoop5_resources', 1335 '%(proc_prefix)sColorNumFactors.dat'%matrix_element.rep_dict),'w') 1336 for ColorLine in CMNum: 1337 CMWriter.write(' '.join(['%d'%C for C in ColorLine])+'\n') 1338 CMWriter.close() 1339 CMWriter=open(pjoin('..','MadLoop5_resources', 1340 '%(proc_prefix)sColorDenomFactors.dat'%matrix_element.rep_dict),'w') 1341 for ColorLine in CMDenom: 1342 CMWriter.write(' '.join(['%d'%C for C in ColorLine])+'\n') 1343 CMWriter.close() 1344 1345 # Write out the helicity configurations 1346 HelConfigs=matrix_element.get_helicity_matrix() 1347 HelConfigWriter=open(pjoin('..','MadLoop5_resources', 1348 '%(proc_prefix)sHelConfigs.dat'%matrix_element.rep_dict),'w') 1349 for HelConfig in HelConfigs: 1350 HelConfigWriter.write(' '.join(['%d'%H for H in HelConfig])+'\n') 1351 HelConfigWriter.close() 1352 1353 # Extract helas calls 1354 loop_amp_helas_calls = fortran_model.get_loop_amp_helas_calls(\ 1355 matrix_element) 1356 # The proc_prefix must be replaced 1357 loop_amp_helas_calls = [lc % matrix_element.rep_dict 1358 for lc in loop_amp_helas_calls] 1359 1360 born_ct_helas_calls, UVCT_helas_calls = \ 1361 fortran_model.get_born_ct_helas_calls(matrix_element) 1362 # In the default output, we do not need to separate these two kind of 1363 # contributions 1364 born_ct_helas_calls = born_ct_helas_calls + UVCT_helas_calls 1365 file = open(os.path.join(self.template_dir,\ 1366 1367 'loop_matrix_standalone.inc')).read() 1368 1369 if matrix_element.get('processes')[0].get('has_born'): 1370 toBeRepaced='loop_helas_calls' 1371 else: 1372 toBeRepaced='loop_induced_helas_calls' 1373 1374 # Decide here wether we need to split the loop_matrix.f file or not. 1375 if (not noSplit and (len(matrix_element.get_all_amplitudes())>1000)): 1376 file=self.split_HELASCALLS(writer,replace_dict,\ 1377 'helas_calls_split.inc',file,born_ct_helas_calls,\ 1378 'born_ct_helas_calls','helas_calls_ampb') 1379 file=self.split_HELASCALLS(writer,replace_dict,\ 1380 'helas_calls_split.inc',file,loop_amp_helas_calls,\ 1381 toBeRepaced,'helas_calls_ampl') 1382 else: 1383 replace_dict['born_ct_helas_calls']='\n'.join(born_ct_helas_calls) 1384 replace_dict[toBeRepaced]='\n'.join(loop_amp_helas_calls) 1385 1386 file = file % replace_dict 1387 1388 loop_calls_finder = re.compile(r'^\s*CALL\S*LOOP\S*') 1389 n_loop_calls = len(filter(lambda call: 1390 not loop_calls_finder.match(call) is None, loop_amp_helas_calls)) 1391 if writer: 1392 # Write the file 1393 writer.writelines(file) 1394 return n_loop_calls 1395 else: 1396 # Return it to be written along with the others 1397 return n_loop_calls, file
1398
1399 - def write_bornmatrix(self, writer, matrix_element, fortran_model):
1400 """Create the born_matrix.f file for the born process as for a standard 1401 tree-level computation.""" 1402 1403 if not matrix_element.get('processes') or \ 1404 not matrix_element.get('diagrams'): 1405 return 0 1406 1407 if not isinstance(writer, writers.FortranWriter): 1408 raise writers.FortranWriter.FortranWriterError(\ 1409 "writer not FortranWriter") 1410 1411 # For now, we can use the exact same treatment as for tree-level 1412 # computations by redefining here a regular HelasMatrixElementf or the 1413 # born process. 1414 # It is important to make a deepcopy, as we don't want any possible 1415 # treatment on the objects of the bornME to have border effects on 1416 # the content of the LoopHelasMatrixElement object. 1417 bornME = helas_objects.HelasMatrixElement() 1418 for prop in bornME.keys(): 1419 bornME.set(prop,copy.deepcopy(matrix_element.get(prop))) 1420 bornME.set('base_amplitude',None,force=True) 1421 bornME.set('diagrams',copy.deepcopy(\ 1422 matrix_element.get_born_diagrams())) 1423 bornME.set('color_basis',copy.deepcopy(\ 1424 matrix_element.get('born_color_basis'))) 1425 bornME.set('color_matrix',copy.deepcopy(\ 1426 color_amp.ColorMatrix(bornME.get('color_basis')))) 1427 # This is to decide wether once to reuse old wavefunction to store new 1428 # ones (provided they are not used further in the code.) 1429 bornME.optimization = True 1430 return super(LoopProcessExporterFortranSA,self).write_matrix_element_v4( 1431 writer, bornME, fortran_model, 1432 proc_prefix=matrix_element.rep_dict['proc_prefix'])
1433
1434 - def write_born_amps_and_wfs(self, writer, matrix_element, fortran_model,\ 1435 noSplit=False):
1436 """ Writes out the code for the subroutine MP_BORN_AMPS_AND_WFS which 1437 computes just the external wavefunction and born amplitudes in 1438 multiple precision. """ 1439 1440 if not matrix_element.get('processes') or \ 1441 not matrix_element.get('diagrams'): 1442 return 0 1443 1444 replace_dict = copy.copy(matrix_element.rep_dict) 1445 1446 # For the wavefunction copy, check what suffix is needed for the W array 1447 if matrix_element.get('processes')[0].get('has_born'): 1448 replace_dict['h_w_suffix']=',H' 1449 else: 1450 replace_dict['h_w_suffix']='' 1451 1452 # Extract helas calls 1453 born_amps_and_wfs_calls , uvct_amp_calls = \ 1454 fortran_model.get_born_ct_helas_calls(matrix_element, include_CT=True) 1455 # In the default output, these two kind of contributions do not need to 1456 # be differentiated 1457 born_amps_and_wfs_calls = born_amps_and_wfs_calls + uvct_amp_calls 1458 1459 # Turn these HELAS calls to the multiple-precision version of the HELAS 1460 # subroutines. 1461 self.turn_to_mp_calls(born_amps_and_wfs_calls) 1462 1463 file = open(os.path.join(self.template_dir,\ 1464 'mp_born_amps_and_wfs.inc')).read() 1465 # Decide here wether we need to split the loop_matrix.f file or not. 1466 if (not noSplit and (len(matrix_element.get_all_amplitudes())>2000)): 1467 file=self.split_HELASCALLS(writer,replace_dict,\ 1468 'mp_helas_calls_split.inc',file,\ 1469 born_amps_and_wfs_calls,'born_amps_and_wfs_calls',\ 1470 'mp_helas_calls') 1471 else: 1472 replace_dict['born_amps_and_wfs_calls']=\ 1473 '\n'.join(born_amps_and_wfs_calls) 1474 1475 file = file % replace_dict 1476 if writer: 1477 # Write the file 1478 writer.writelines(file) 1479 else: 1480 # Return it to be written along with the others 1481 return file 1482 1483 #=============================================================================== 1484 # LoopProcessOptimizedExporterFortranSA 1485 #=============================================================================== 1486
1487 -class LoopProcessOptimizedExporterFortranSA(LoopProcessExporterFortranSA):
1488 """Class to take care of exporting a set of loop matrix elements in the 1489 Fortran format which exploits the Pozzorini method of representing 1490 the loop numerators as polynomial to render its evaluations faster.""" 1491 1492 template_dir=os.path.join(_file_path,'iolibs/template_files/loop_optimized') 1493 # The option below controls wether one wants to group together in one single 1494 # CutTools/TIR call the loops with same denominator structure 1495 forbid_loop_grouping = False 1496 1497 # List of potential TIR library one wants to link to 1498 all_tir=['pjfry','iregi','golem'] 1499
1500 - def __init__(self, mgme_dir="", dir_path = "", opt=None):
1501 """Initiate the LoopProcessOptimizedExporterFortranSA with directory 1502 information on where to find all the loop-related source files, 1503 like CutTools and TIR""" 1504 1505 super(LoopProcessOptimizedExporterFortranSA,self).__init__(mgme_dir, 1506 dir_path, opt) 1507 1508 # TIR available ones 1509 self.tir_available_dict={'pjfry':True,'iregi':True,'golem':True} 1510 1511 for tir in self.all_tir: 1512 tir_dir="%s_dir"%tir 1513 if tir_dir in self.opt: 1514 setattr(self,tir_dir,self.opt[tir_dir]) 1515 else: 1516 setattr(self,tir_dir,'')
1517
1518 - def copy_v4template(self, modelname = ''):
1519 """Additional actions needed to setup the Template. 1520 """ 1521 1522 super(LoopProcessOptimizedExporterFortranSA, self).copy_v4template( 1523 modelname) 1524 1525 self.loop_optimized_additional_template_setup()
1526
1527 - def get_context(self,matrix_element, **opts):
1528 """ Additional contextual information which needs to be created for 1529 the optimized output.""" 1530 1531 context = LoopProcessExporterFortranSA.get_context(self, matrix_element, 1532 **opts) 1533 1534 for tir in self.all_tir: 1535 context['%s_available'%tir]=self.tir_available_dict[tir] 1536 # safety check 1537 if tir not in ['golem','pjfry','iregi']: 1538 raise MadGraph5Error,"%s was not a TIR currently interfected."%tir_name 1539 1540 return context
1541
1543 """ Perform additional actions specific for this class when setting 1544 up the template with the copy_v4template function.""" 1545 1546 # We must link the TIR to the Library folder of the active Template 1547 link_tir_libs=[] 1548 tir_libs=[] 1549 tir_include=[] 1550 # special for PJFry++ 1551 link_pjfry_lib="" 1552 pjfry_lib="" 1553 pjdir="" 1554 for tir in self.all_tir: 1555 tir_dir="%s_dir"%tir 1556 libpath=getattr(self,tir_dir) 1557 libname="lib%s.a"%tir 1558 tir_name=tir 1559 libpath = self.link_TIR(os.path.join(self.dir_path, 'lib'), 1560 libpath,libname,tir_name=tir_name) 1561 setattr(self,tir_dir,libpath) 1562 if libpath != "": 1563 if tir in ['pjfry','golem']: 1564 # Apparently it is necessary to link against the original 1565 # location of the pjfry library, so it needs a special treatment. 1566 link_tir_libs.append('-L%s/ -l%s'%(libpath,tir)) 1567 tir_libs.append('%s/lib%s.$(libext)'%(libpath,tir)) 1568 if tir=='golem': 1569 trgt_path = pjoin(os.path.dirname(libpath),'include') 1570 golem_include = misc.find_includes_path(trgt_path,'.mod') 1571 if golem_include is None: 1572 logger.error( 1573 'Could not find the include directory for golem, looking in %s.\n' % str(trgt_path)+ 1574 'Generation carries on but you will need to edit the include path by hand in the makefiles.') 1575 golem_include = '<Not_found_define_it_yourself>' 1576 tir_include.append('-I %s'%str(golem_include)) 1577 # To be able to easily compile a MadLoop library using 1578 # makefiles built outside of the MG5_aMC framework 1579 # (such as what is done with the Sherpa interface), we 1580 # place here an easy handle on the golem includes 1581 ln(golem_include, starting_dir=pjoin(self.dir_path,'lib'), 1582 name='golem95_include',abspath=True) 1583 ln(libpath, starting_dir=pjoin(self.dir_path,'lib'), 1584 name='golem95_lib',abspath=True) 1585 1586 else : 1587 link_tir_libs.append('-l%s'%tir) 1588 tir_libs.append('$(LIBDIR)lib%s.$(libext)'%tir) 1589 1590 MadLoop_makefile_definitions = pjoin(self.dir_path,'SubProcesses', 1591 'MadLoop_makefile_definitions') 1592 if os.path.isfile(MadLoop_makefile_definitions): 1593 os.remove(MadLoop_makefile_definitions) 1594 1595 calls = self.write_loop_makefile_definitions( 1596 writers.MakefileWriter(MadLoop_makefile_definitions), 1597 link_tir_libs,tir_libs, tir_include=tir_include)
1598 1609 1610 1721
1722 - def set_group_loops(self, matrix_element):
1723 """ Decides whether we must group loops or not for this matrix element""" 1724 1725 # Decide if loops sharing same denominator structures have to be grouped 1726 # together or not. 1727 if self.forbid_loop_grouping: 1728 self.group_loops = False 1729 else: 1730 self.group_loops = (not self.get_context(matrix_element)['ComputeColorFlows'])\ 1731 and matrix_element.get('processes')[0].get('has_born') 1732 1733 return self.group_loops
1734
1735 - def write_loop_matrix_element_v4(self, writer, matrix_element, fortran_model, 1736 group_number = None, proc_id = None, config_map = None):
1737 """ Writes loop_matrix.f, CT_interface.f,TIR_interface.f,GOLEM_inteface.f 1738 and loop_num.f only but with the optimized FortranModel. 1739 The arguments group_number and proc_id are just for the LoopInduced 1740 output with MadEvent and only used in get_ME_identifier.""" 1741 1742 # Warn the user that the 'matrix' output where all relevant code is 1743 # put together in a single file is not supported in this loop output. 1744 if writer: 1745 raise MadGraph5Error, 'Matrix output mode no longer supported.' 1746 1747 if not isinstance(fortran_model,\ 1748 helas_call_writers.FortranUFOHelasCallWriter): 1749 raise MadGraph5Error, 'The optimized loop fortran output can only'+\ 1750 ' work with a UFO Fortran model' 1751 OptimizedFortranModel=\ 1752 helas_call_writers.FortranUFOHelasCallWriterOptimized(\ 1753 fortran_model.get('model'),False) 1754 1755 1756 if not matrix_element.get('processes')[0].get('has_born') and \ 1757 not self.compute_color_flows: 1758 logger.debug("Color flows will be employed despite the option"+\ 1759 " 'loop_color_flows' being set to False because it is necessary"+\ 1760 " for optimizations.") 1761 1762 # Compute the analytical information of the loop wavefunctions in the 1763 # loop helas matrix elements using the cached aloha model to reuse 1764 # as much as possible the aloha computations already performed for 1765 # writing out the aloha fortran subroutines. 1766 matrix_element.compute_all_analytic_information( 1767 self.get_aloha_model(matrix_element.get('processes')[0].get('model'))) 1768 1769 self.set_group_loops(matrix_element) 1770 1771 # Initialize a general replacement dictionary with entries common to 1772 # many files generated here. 1773 matrix_element.rep_dict = LoopProcessExporterFortranSA.\ 1774 generate_general_replace_dict(self, matrix_element, 1775 group_number = group_number, proc_id = proc_id) 1776 1777 # and those specific to the optimized output 1778 self.set_optimized_output_specific_replace_dict_entries(matrix_element) 1779 1780 # Create the necessary files for the loop matrix element subroutine 1781 proc_prefix_writer = writers.FortranWriter('proc_prefix.txt','w') 1782 proc_prefix_writer.write(matrix_element.rep_dict['proc_prefix']) 1783 proc_prefix_writer.close() 1784 1785 filename = 'loop_matrix.f' 1786 calls = self.write_loopmatrix(writers.FortranWriter(filename), 1787 matrix_element, 1788 OptimizedFortranModel) 1789 1790 filename = 'check_sa.f' 1791 self.write_check_sa(writers.FortranWriter(filename),matrix_element) 1792 1793 filename = 'polynomial.f' 1794 calls = self.write_polynomial_subroutines( 1795 writers.FortranWriter(filename), 1796 matrix_element) 1797 1798 filename = 'improve_ps.f' 1799 calls = self.write_improve_ps(writers.FortranWriter(filename), 1800 matrix_element) 1801 1802 filename = 'CT_interface.f' 1803 self.write_CT_interface(writers.FortranWriter(filename),\ 1804 matrix_element) 1805 1806 filename = 'TIR_interface.f' 1807 self.write_TIR_interface(writers.FortranWriter(filename), 1808 matrix_element) 1809 1810 if 'golem' in self.tir_available_dict and self.tir_available_dict['golem']: 1811 filename = 'GOLEM_interface.f' 1812 self.write_GOLEM_interface(writers.FortranWriter(filename), 1813 matrix_element) 1814 1815 filename = 'loop_num.f' 1816 self.write_loop_num(writers.FortranWriter(filename),\ 1817 matrix_element,OptimizedFortranModel) 1818 1819 filename = 'mp_compute_loop_coefs.f' 1820 self.write_mp_compute_loop_coefs(writers.FortranWriter(filename),\ 1821 matrix_element,OptimizedFortranModel) 1822 1823 if self.get_context(matrix_element)['ComputeColorFlows']: 1824 filename = 'compute_color_flows.f' 1825 self.write_compute_color_flows(writers.FortranWriter(filename), 1826 matrix_element, config_map = config_map) 1827 1828 # Extract number of external particles 1829 (nexternal, ninitial) = matrix_element.get_nexternal_ninitial() 1830 filename = 'nexternal.inc' 1831 self.write_nexternal_file(writers.FortranWriter(filename), 1832 nexternal, ninitial) 1833 1834 if self.get_context(matrix_element)['TIRCaching']: 1835 filename = 'tir_cache_size.inc' 1836 self.write_tir_cache_size_include(writers.FortranWriter(filename)) 1837 1838 return calls
1839
1840 - def set_optimized_output_specific_replace_dict_entries(self, matrix_element):
1841 """ Specify the entries of the replacement dictionary which are specific 1842 to the optimized output and only relevant to it (the more general entries 1843 are set in the the mother class LoopProcessExporterFortranSA.""" 1844 1845 max_loop_rank=matrix_element.get_max_loop_rank() 1846 matrix_element.rep_dict['maxrank']=max_loop_rank 1847 matrix_element.rep_dict['loop_max_coefs']=\ 1848 q_polynomial.get_number_of_coefs_for_rank(max_loop_rank) 1849 max_loop_vertex_rank=matrix_element.get_max_loop_vertex_rank() 1850 matrix_element.rep_dict['vertex_max_coefs']=\ 1851 q_polynomial.get_number_of_coefs_for_rank(max_loop_vertex_rank) 1852 1853 matrix_element.rep_dict['nloopwavefuncs']=\ 1854 matrix_element.get_number_of_loop_wavefunctions() 1855 max_spin=matrix_element.get_max_loop_particle_spin() 1856 if max_spin>3: 1857 raise MadGraph5Error, "ML5 can only handle loop particles with"+\ 1858 " spin 1 at most" 1859 matrix_element.rep_dict['max_lwf_size']=4 1860 matrix_element.rep_dict['nloops']=len(\ 1861 [1 for ldiag in matrix_element.get_loop_diagrams() for \ 1862 lamp in ldiag.get_loop_amplitudes()]) 1863 1864 if self.set_group_loops(matrix_element): 1865 matrix_element.rep_dict['nloop_groups']=\ 1866 len(matrix_element.get('loop_groups')) 1867 else: 1868 matrix_element.rep_dict['nloop_groups']=\ 1869 matrix_element.rep_dict['nloops']
1870
1871 - def write_loop_num(self, writer, matrix_element,fortran_model):
1872 """ Create the file containing the core subroutine called by CutTools 1873 which contains the Helas calls building the loop""" 1874 1875 replace_dict=copy.copy(matrix_element.rep_dict) 1876 1877 file = open(os.path.join(self.template_dir,'loop_num.inc')).read() 1878 file = file % replace_dict 1879 writer.writelines(file,context=self.get_context(matrix_element))
1880
1881 - def write_CT_interface(self, writer, matrix_element):
1882 """ We can re-use the mother one for the loop optimized output.""" 1883 LoopProcessExporterFortranSA.write_CT_interface(\ 1884 self, writer, matrix_element,optimized_output=True)
1885
1886 - def write_TIR_interface(self, writer, matrix_element):
1887 """ Create the file TIR_interface.f which does NOT contain the subroutine 1888 defining the loop HELAS-like calls along with the general interfacing 1889 subroutine. """ 1890 1891 # First write TIR_interface which interfaces MG5 with TIR. 1892 replace_dict=copy.copy(matrix_element.rep_dict) 1893 1894 file = open(os.path.join(self.template_dir,'TIR_interface.inc')).read() 1895 1896 file = file % replace_dict 1897 1898 FPR = q_polynomial.FortranPolynomialRoutines( 1899 replace_dict['maxrank'],coef_format=replace_dict['complex_dp_format'],\ 1900 sub_prefix=replace_dict['proc_prefix']) 1901 if self.tir_available_dict['pjfry']: 1902 file += '\n\n'+FPR.write_pjfry_mapping() 1903 if self.tir_available_dict['iregi']: 1904 file += '\n\n'+FPR.write_iregi_mapping() 1905 1906 if writer: 1907 writer.writelines(file,context=self.get_context(matrix_element)) 1908 else: 1909 return file
1910
1911 - def write_GOLEM_interface(self, writer, matrix_element):
1912 """ Create the file GOLEM_interface.f which does NOT contain the subroutine 1913 defining the loop HELAS-like calls along with the general interfacing 1914 subroutine. """ 1915 1916 # First write GOLEM_interface which interfaces MG5 with TIR. 1917 replace_dict=copy.copy(matrix_element.rep_dict) 1918 1919 # We finalize TIR result differently wether we used the built-in 1920 # squaring against the born. 1921 if not self.get_context(matrix_element)['AmplitudeReduction']: 1922 replace_dict['loop_induced_sqsoindex']=',SQSOINDEX' 1923 else: 1924 replace_dict['loop_induced_sqsoindex']='' 1925 1926 file = open(os.path.join(self.template_dir,'GOLEM_interface.inc')).read() 1927 1928 file = file % replace_dict 1929 1930 FPR = q_polynomial.FortranPolynomialRoutines(replace_dict['maxrank'],\ 1931 coef_format=replace_dict['complex_dp_format'],\ 1932 sub_prefix=replace_dict['proc_prefix']) 1933 1934 file += '\n\n'+FPR.write_golem95_mapping() 1935 1936 if writer: 1937 writer.writelines(file,context=self.get_context(matrix_element)) 1938 else: 1939 return file
1940
1941 - def write_polynomial_subroutines(self,writer,matrix_element):
1942 """ Subroutine to create all the subroutines relevant for handling 1943 the polynomials representing the loop numerator """ 1944 1945 # First create 'coef_specs.inc' 1946 IncWriter=writers.FortranWriter('coef_specs.inc','w') 1947 IncWriter.writelines("""INTEGER MAXLWFSIZE 1948 PARAMETER (MAXLWFSIZE=%(max_lwf_size)d) 1949 INTEGER LOOP_MAXCOEFS 1950 PARAMETER (LOOP_MAXCOEFS=%(loop_max_coefs)d) 1951 INTEGER VERTEXMAXCOEFS 1952 PARAMETER (VERTEXMAXCOEFS=%(vertex_max_coefs)d)"""\ 1953 %matrix_element.rep_dict) 1954 IncWriter.close() 1955 1956 # List of all subroutines to place there 1957 subroutines=[] 1958 1959 # Start from the routine in the template 1960 replace_dict = copy.copy(matrix_element.rep_dict) 1961 dp_routine = open(os.path.join(self.template_dir,'polynomial.inc')).read() 1962 mp_routine = open(os.path.join(self.template_dir,'polynomial.inc')).read() 1963 # The double precision version of the basic polynomial routines, such as 1964 # create_loop_coefs 1965 replace_dict['complex_format'] = replace_dict['complex_dp_format'] 1966 replace_dict['real_format'] = replace_dict['real_dp_format'] 1967 replace_dict['mp_prefix'] = '' 1968 replace_dict['kind'] = 8 1969 replace_dict['zero_def'] = '0.0d0' 1970 replace_dict['one_def'] = '1.0d0' 1971 dp_routine = dp_routine % replace_dict 1972 # The quadruple precision version of the basic polynomial routines 1973 replace_dict['complex_format'] = replace_dict['complex_mp_format'] 1974 replace_dict['real_format'] = replace_dict['real_mp_format'] 1975 replace_dict['mp_prefix'] = 'MP_' 1976 replace_dict['kind'] = 16 1977 replace_dict['zero_def'] = '0.0e0_16' 1978 replace_dict['one_def'] = '1.0e0_16' 1979 mp_routine = mp_routine % replace_dict 1980 subroutines.append(dp_routine) 1981 subroutines.append(mp_routine) 1982 1983 # Initialize the polynomial routine writer 1984 poly_writer=q_polynomial.FortranPolynomialRoutines( 1985 matrix_element.get_max_loop_rank(), 1986 sub_prefix=replace_dict['proc_prefix']) 1987 mp_poly_writer=q_polynomial.FortranPolynomialRoutines( 1988 matrix_element.get_max_loop_rank(),coef_format='complex*32', 1989 sub_prefix='MP_'+replace_dict['proc_prefix']) 1990 # The eval subroutine 1991 subroutines.append(poly_writer.write_polynomial_evaluator()) 1992 subroutines.append(mp_poly_writer.write_polynomial_evaluator()) 1993 # The add coefs subroutine 1994 subroutines.append(poly_writer.write_add_coefs()) 1995 subroutines.append(mp_poly_writer.write_add_coefs()) 1996 # The merging one for creating the loop coefficients 1997 subroutines.append(poly_writer.write_wl_merger()) 1998 subroutines.append(mp_poly_writer.write_wl_merger()) 1999 # Now the udpate subroutines 2000 for wl_update in matrix_element.get_used_wl_updates(): 2001 subroutines.append(poly_writer.write_wl_updater(\ 2002 wl_update[0],wl_update[1])) 2003 subroutines.append(mp_poly_writer.write_wl_updater(\ 2004 wl_update[0],wl_update[1])) 2005 writer.writelines('\n\n'.join(subroutines), 2006 context=self.get_context(matrix_element))
2007
2008 - def write_mp_compute_loop_coefs(self, writer, matrix_element, fortran_model, \ 2009 noSplit=False):
2010 """Create the write_mp_compute_loop_coefs.f file.""" 2011 2012 if not matrix_element.get('processes') or \ 2013 not matrix_element.get('diagrams'): 2014 return 0 2015 2016 # Set lowercase/uppercase Fortran code 2017 2018 writers.FortranWriter.downcase = False 2019 2020 replace_dict = copy.copy(matrix_element.rep_dict) 2021 2022 # Extract helas calls 2023 squared_orders = matrix_element.get_squared_order_contribs() 2024 split_orders = matrix_element.get('processes')[0].get('split_orders') 2025 2026 born_ct_helas_calls , uvct_helas_calls = \ 2027 fortran_model.get_born_ct_helas_calls(matrix_element, 2028 squared_orders=squared_orders, split_orders=split_orders) 2029 self.turn_to_mp_calls(born_ct_helas_calls) 2030 self.turn_to_mp_calls(uvct_helas_calls) 2031 coef_construction, coef_merging = fortran_model.get_coef_construction_calls(\ 2032 matrix_element,group_loops=self.group_loops, 2033 squared_orders=squared_orders,split_orders=split_orders) 2034 # The proc_prefix must be replaced 2035 coef_construction = [c % matrix_element.rep_dict for c 2036 in coef_construction] 2037 self.turn_to_mp_calls(coef_construction) 2038 self.turn_to_mp_calls(coef_merging) 2039 2040 file = open(os.path.join(self.template_dir,\ 2041 'mp_compute_loop_coefs.inc')).read() 2042 2043 # Setup the contextual environment which is used in the splitting 2044 # functions below 2045 context = self.get_context(matrix_element) 2046 # Decide here wether we need to split the loop_matrix.f file or not. 2047 # 200 is reasonable but feel free to change it. 2048 if (not noSplit and (len(matrix_element.get_all_amplitudes())>200)): 2049 file=self.split_HELASCALLS(writer,replace_dict,\ 2050 'mp_helas_calls_split.inc',file,born_ct_helas_calls,\ 2051 'mp_born_ct_helas_calls','mp_helas_calls_ampb', 2052 required_so_broadcaster = 'MP_CT_REQ_SO_DONE', 2053 continue_label = 2000,context=context) 2054 file=self.split_HELASCALLS(writer,replace_dict,\ 2055 'mp_helas_calls_split.inc',file,uvct_helas_calls,\ 2056 'mp_uvct_helas_calls','mp_helas_calls_uvct', 2057 required_so_broadcaster = 'MP_UVCT_REQ_SO_DONE', 2058 continue_label = 3000,context=context) 2059 file=self.split_HELASCALLS(writer,replace_dict,\ 2060 'mp_helas_calls_split.inc',file,coef_construction,\ 2061 'mp_coef_construction','mp_coef_construction', 2062 required_so_broadcaster = 'MP_LOOP_REQ_SO_DONE', 2063 continue_label = 4000,context=context) 2064 else: 2065 replace_dict['mp_born_ct_helas_calls']='\n'.join(born_ct_helas_calls) 2066 replace_dict['mp_uvct_helas_calls']='\n'.join(uvct_helas_calls) 2067 replace_dict['mp_coef_construction']='\n'.join(coef_construction) 2068 2069 replace_dict['mp_coef_merging']='\n'.join(coef_merging) 2070 2071 file = file % replace_dict 2072 2073 # Write the file 2074 writer.writelines(file,context=context)
2075
2076 - def write_color_matrix_data_file(self, writer, col_matrix):
2077 """Writes out the files (Loop|Born)ColorFlowMatrix.dat corresponding 2078 to the color coefficients for JAMP(L|B)*JAMP(L|B).""" 2079 2080 res = [] 2081 for line in range(len(col_matrix._col_basis1)): 2082 numerators = [] 2083 denominators = [] 2084 for row in range(len(col_matrix._col_basis2)): 2085 coeff = col_matrix.col_matrix_fixed_Nc[(line,row)] 2086 numerators.append('%6r'%coeff[0].numerator) 2087 denominators.append('%6r'%( 2088 coeff[0].denominator*(-1 if coeff[1] else 1))) 2089 res.append(' '.join(numerators)) 2090 res.append(' '.join(denominators)) 2091 2092 res.append('EOF') 2093 2094 writer.writelines('\n'.join(res))
2095
2096 - def write_color_flow_coefs_data_file(self, writer, color_amplitudes, 2097 color_basis):
2098 """ Writes the file '(Loop|Born)ColorFlowCoefs.dat using the coefficients 2099 list of the color_amplitudes in the argument of this function.""" 2100 2101 my_cs = color.ColorString() 2102 2103 res = [] 2104 2105 for jamp_number, coeff_list in enumerate(color_amplitudes): 2106 my_cs.from_immutable(sorted(color_basis.keys())[jamp_number]) 2107 # Order the ColorString so that its ordering is canonical. 2108 ordered_cs = color.ColorFactor([my_cs]).full_simplify()[0] 2109 res.append('%d # Coefficient for flow number %d with expr. %s'\ 2110 %(len(coeff_list), jamp_number+1, repr(ordered_cs))) 2111 # A line element is a tuple (numerator, denominator, amplitude_id) 2112 line_element = [] 2113 2114 for (coefficient, amp_number) in coeff_list: 2115 coef = self.cat_coeff(\ 2116 coefficient[0],coefficient[1],coefficient[2],coefficient[3]) 2117 line_element.append((coef[0].numerator, 2118 coef[0].denominator*(-1 if coef[1] else 1),amp_number)) 2119 # Sort them by growing amplitude number 2120 line_element.sort(key=lambda el:el[2]) 2121 2122 for i in range(3): 2123 res.append(' '.join('%6r'%elem[i] for elem in line_element)) 2124 2125 res.append('EOF') 2126 writer.writelines('\n'.join(res))
2127
2128 - def write_compute_color_flows(self, writer, matrix_element, config_map):
2129 """Writes the file compute_color_flows.f which uses the AMPL results 2130 from a common block to project them onto the color flow space so as 2131 to compute the JAMP quantities. For loop induced processes, this file 2132 will also contain a subroutine computing AMPL**2 for madevent 2133 multichanneling.""" 2134 2135 loop_col_amps = matrix_element.get_loop_color_amplitudes() 2136 matrix_element.rep_dict['nLoopFlows'] = len(loop_col_amps) 2137 2138 dat_writer = open(pjoin('..','MadLoop5_resources', 2139 '%(proc_prefix)sLoopColorFlowCoefs.dat' 2140 %matrix_element.rep_dict),'w') 2141 self.write_color_flow_coefs_data_file(dat_writer, 2142 loop_col_amps, matrix_element.get('loop_color_basis')) 2143 dat_writer.close() 2144 2145 dat_writer = open(pjoin('..','MadLoop5_resources', 2146 '%(proc_prefix)sLoopColorFlowMatrix.dat' 2147 %matrix_element.rep_dict),'w') 2148 self.write_color_matrix_data_file(dat_writer, 2149 matrix_element.get('color_matrix')) 2150 dat_writer.close() 2151 2152 if matrix_element.get('processes')[0].get('has_born'): 2153 born_col_amps = matrix_element.get_born_color_amplitudes() 2154 matrix_element.rep_dict['nBornFlows'] = len(born_col_amps) 2155 dat_writer = open(pjoin('..','MadLoop5_resources', 2156 '%(proc_prefix)sBornColorFlowCoefs.dat' 2157 %matrix_element.rep_dict),'w') 2158 self.write_color_flow_coefs_data_file(dat_writer, 2159 born_col_amps, matrix_element.get('loop_color_basis')) 2160 dat_writer.close() 2161 2162 dat_writer = open(pjoin('..','MadLoop5_resources', 2163 '%(proc_prefix)sBornColorFlowMatrix.dat' 2164 %matrix_element.rep_dict),'w') 2165 self.write_color_matrix_data_file(dat_writer, 2166 color_amp.ColorMatrix(matrix_element.get('born_color_basis'))) 2167 dat_writer.close() 2168 else: 2169 matrix_element.rep_dict['nBornFlows'] = 0 2170 2171 replace_dict = copy.copy(matrix_element.rep_dict) 2172 2173 # The following variables only have to be defined for the LoopInduced 2174 # output for madevent. 2175 if self.get_context(matrix_element)['MadEventOutput']: 2176 self.get_amp2_lines(matrix_element, replace_dict, config_map) 2177 else: 2178 replace_dict['config_map_definition'] = '' 2179 replace_dict['config_index_map_definition'] = '' 2180 replace_dict['nmultichannels'] = 0 2181 replace_dict['nmultichannel_configs'] = 0 2182 2183 # The nmultichannels entry will be used in the matrix<i> wrappers as 2184 # well, so we add it to the general_replace_dict too. 2185 matrix_element.rep_dict['nmultichannels'] = \ 2186 replace_dict['nmultichannels'] 2187 matrix_element.rep_dict['nmultichannel_configs'] = \ 2188 replace_dict['nmultichannel_configs'] 2189 2190 2191 file = open(os.path.join(self.template_dir,\ 2192 'compute_color_flows.inc')).read()%replace_dict 2193 2194 writer.writelines(file,context=self.get_context(matrix_element))
2195
2196 - def fix_coef_specs(self, overall_max_lwf_size, overall_max_loop_vert_rank):
2197 """ If processes with different maximum loop wavefunction size or 2198 different maximum loop vertex rank have to be output together, then 2199 the file 'coef.inc' in the HELAS Source folder must contain the overall 2200 maximum of these quantities. It is not safe though, and the user has 2201 been appropriatly warned at the output stage """ 2202 2203 # Remove the existing link 2204 coef_specs_path=os.path.join(self.dir_path,'Source','DHELAS',\ 2205 'coef_specs.inc') 2206 os.remove(coef_specs_path) 2207 2208 # Replace it by the appropriate value 2209 IncWriter=writers.FortranWriter(coef_specs_path,'w') 2210 IncWriter.writelines("""INTEGER MAXLWFSIZE 2211 PARAMETER (MAXLWFSIZE=%(max_lwf_size)d) 2212 INTEGER VERTEXMAXCOEFS 2213 PARAMETER (VERTEXMAXCOEFS=%(vertex_max_coefs)d)"""\ 2214 %{'max_lwf_size':overall_max_lwf_size, 2215 'vertex_max_coefs':overall_max_loop_vert_rank}) 2216 IncWriter.close()
2217
2218 - def setup_check_sa_replacement_dictionary(self, matrix_element, \ 2219 split_orders,squared_orders,amps_orders):
2220 """ Sets up the replacement dictionary for the writeout of the steering 2221 file check_sa.f""" 2222 if len(squared_orders)<1: 2223 matrix_element.rep_dict['print_so_loop_results']=\ 2224 "write(*,*) 'No split orders defined.'" 2225 elif len(squared_orders)==1: 2226 matrix_element.rep_dict['set_coupling_target']='' 2227 matrix_element.rep_dict['print_so_loop_results']=\ 2228 "write(*,*) 'All loop contributions are of split orders (%s)'"%( 2229 ' '.join(['%s=%d'%(split_orders[i],squared_orders[0][i]) \ 2230 for i in range(len(split_orders))])) 2231 else: 2232 matrix_element.rep_dict['set_coupling_target']='\n'.join([ 2233 '# Here we leave the default target squared split order to -1, meaning that we'+ 2234 ' aim at computing all individual contributions. You can choose otherwise.', 2235 'call %(proc_prefix)sSET_COUPLINGORDERS_TARGET(-1)'%matrix_element.rep_dict]) 2236 matrix_element.rep_dict['print_so_loop_results'] = '\n'.join([ 2237 '\n'.join(["write(*,*) '%dL) Loop ME for orders (%s) :'"%((j+1),(' '.join( 2238 ['%s=%d'%(split_orders[i],so[i]) for i in range(len(split_orders))]))), 2239 "IF (PREC_FOUND(%d).NE.-1.0d0) THEN"%(j+1), 2240 "write(*,*) ' > accuracy = ',PREC_FOUND(%d)"%(j+1), 2241 "ELSE", 2242 "write(*,*) ' > accuracy = NA'", 2243 "ENDIF", 2244 "write(*,*) ' > finite = ',MATELEM(1,%d)"%(j+1), 2245 "write(*,*) ' > 1eps = ',MATELEM(2,%d)"%(j+1), 2246 "write(*,*) ' > 2eps = ',MATELEM(3,%d)"%(j+1) 2247 ]) for j, so in enumerate(squared_orders)]) 2248 matrix_element.rep_dict['write_so_loop_results'] = '\n'.join( 2249 ["write (69,*) 'Split_Orders_Names %s'"%(' '.join(split_orders))]+ 2250 ['\n'.join([ 2251 "write (69,*) 'Loop_SO_Results %s'"%(' '.join( 2252 ['%d'%so_value for so_value in so])), 2253 "write (69,*) 'SO_Loop ACC ',PREC_FOUND(%d)"%(j+1), 2254 "write (69,*) 'SO_Loop FIN ',MATELEM(1,%d)"%(j+1), 2255 "write (69,*) 'SO_Loop 1EPS ',MATELEM(2,%d)"%(j+1), 2256 "write (69,*) 'SO_Loop 2EPS ',MATELEM(3,%d)"%(j+1), 2257 ]) for j, so in enumerate(squared_orders)]) 2258 2259 # We must reconstruct here the born squared orders. 2260 squared_born_so_orders = [] 2261 for i, amp_order in enumerate(amps_orders['born_amp_orders']): 2262 for j in range(0,i+1): 2263 key = tuple([ord1 + ord2 for ord1,ord2 in \ 2264 zip(amp_order[0],amps_orders['born_amp_orders'][j][0])]) 2265 if not key in squared_born_so_orders: 2266 squared_born_so_orders.append(key) 2267 if len(squared_born_so_orders)<1: 2268 matrix_element.rep_dict['print_so_born_results'] = '' 2269 elif len(squared_born_so_orders)==1: 2270 matrix_element.rep_dict['print_so_born_results'] = \ 2271 "write(*,*) 'All Born contributions are of split orders (%s)'"%( 2272 ' '.join(['%s=%d'%(split_orders[i],squared_born_so_orders[0][i]) 2273 for i in range(len(split_orders))])) 2274 else: 2275 matrix_element.rep_dict['print_so_born_results'] = '\n'.join([ 2276 "write(*,*) '%dB) Born ME for orders (%s) = ',MATELEM(0,%d)"%(j+1,' '.join( 2277 ['%s=%d'%(split_orders[i],so[i]) for i in range(len(split_orders))]),j+1) 2278 for j, so in enumerate(squared_born_so_orders)]) 2279 matrix_element.rep_dict['write_so_born_results'] = '\n'.join( 2280 ['\n'.join([ 2281 "write (69,*) 'Born_SO_Results %s'"%(' '.join( 2282 ['%d'%so_value for so_value in so])), 2283 "write (69,*) 'SO_Born BORN ',MATELEM(0,%d)"%(j+1), 2284 ]) for j, so in enumerate(squared_born_so_orders)]) 2285 2286 # Add a bottom bar to both print_so_[loop|born]_results 2287 matrix_element.rep_dict['print_so_born_results'] += \ 2288 '\nwrite (*,*) "---------------------------------"' 2289 matrix_element.rep_dict['print_so_loop_results'] += \ 2290 '\nwrite (*,*) "---------------------------------"'
2291
2292 - def write_tir_cache_size_include(self, writer):
2293 """Write the file 'tir_cache_size.inc' which sets the size of the TIR 2294 cache the the user wishes to employ and the default value for it. 2295 This can have an impact on MadLoop speed when using stability checks 2296 but also impacts in a non-negligible way MadLoop's memory footprint. 2297 It is therefore important that the user can chose its size.""" 2298 2299 # For the standalone optimized output, a size of one is necessary. 2300 # The MadLoop+MadEvent output sets it to 2 because it can gain further 2301 # speed increase with a TIR cache of size 2 due to the structure of the 2302 # calls to MadLoop there. 2303 tir_cach_size = "parameter(TIR_CACHE_SIZE=1)" 2304 writer.writelines(tir_cach_size)
2305
2306 - def write_loopmatrix(self, writer, matrix_element, fortran_model, \ 2307 noSplit=False, write_auxiliary_files=True,):
2308 """Create the loop_matrix.f file.""" 2309 2310 if not matrix_element.get('processes') or \ 2311 not matrix_element.get('diagrams'): 2312 return 0 2313 2314 # Set lowercase/uppercase Fortran code 2315 writers.FortranWriter.downcase = False 2316 2317 # Starting off with the treatment of the split_orders since some 2318 # of the information extracted there will come into the 2319 # general_replace_dict. Split orders are abbreviated SO in all the 2320 # keys of the replacement dictionaries. 2321 2322 # Take care of the split_orders 2323 squared_orders, amps_orders = matrix_element.get_split_orders_mapping() 2324 # Creating here a temporary list containing only the information of 2325 # what are the different squared split orders contributing 2326 # (i.e. not using max_contrib_amp_number and max_contrib_ref_amp_number) 2327 sqso_contribs = [sqso[0] for sqso in squared_orders] 2328 split_orders = matrix_element.get('processes')[0].get('split_orders') 2329 # The entries set in the function below are only for check_sa written 2330 # out in write_loop__matrix_element_v4 (it is however placed here because the 2331 # split order information is only available here). 2332 self.setup_check_sa_replacement_dictionary(matrix_element, 2333 split_orders,sqso_contribs,amps_orders) 2334 2335 # Now recast the split order basis for the loop, born and counterterm 2336 # amplitude into one single splitorderbasis. 2337 overall_so_basis = list(set( 2338 [born_so[0] for born_so in amps_orders['born_amp_orders']]+ 2339 [born_so[0] for born_so in amps_orders['loop_amp_orders']])) 2340 # We must re-sort it to make sure it follows an increasing WEIGHT order 2341 order_hierarchy = matrix_element.get('processes')[0]\ 2342 .get('model').get('order_hierarchy') 2343 if set(order_hierarchy.keys()).union(set(split_orders))==\ 2344 set(order_hierarchy.keys()): 2345 overall_so_basis.sort(key= lambda so: 2346 sum([order_hierarchy[split_orders[i]]*order_power for \ 2347 i, order_power in enumerate(so)])) 2348 2349 # Those are additional entries used throughout the different files of 2350 # MadLoop5 2351 matrix_element.rep_dict['split_order_str_list'] = str(split_orders) 2352 matrix_element.rep_dict['nSO'] = len(split_orders) 2353 matrix_element.rep_dict['nSquaredSO'] = len(sqso_contribs) 2354 matrix_element.rep_dict['nAmpSO'] = len(overall_so_basis) 2355 2356 writers.FortranWriter('nsquaredSO.inc').writelines( 2357 """INTEGER NSQUAREDSO 2358 PARAMETER (NSQUAREDSO=%d)"""%matrix_element.rep_dict['nSquaredSO']) 2359 2360 replace_dict = copy.copy(matrix_element.rep_dict) 2361 # Build the general array mapping the split orders indices to their 2362 # definition 2363 replace_dict['ampsplitorders'] = '\n'.join(self.get_split_orders_lines(\ 2364 overall_so_basis,'AMPSPLITORDERS')) 2365 replace_dict['SquaredSO'] = '\n'.join(self.get_split_orders_lines(\ 2366 sqso_contribs,'SQPLITORDERS')) 2367 2368 # Specify what are the squared split orders selected by the proc def. 2369 replace_dict['chosen_so_configs'] = self.set_chosen_SO_index( 2370 matrix_element.get('processes')[0],sqso_contribs) 2371 2372 # Now we build the different arrays storing the split_orders ID of each 2373 # amp. 2374 ampSO_list=[-1]*sum(len(el[1]) for el in amps_orders['loop_amp_orders']) 2375 for SO in amps_orders['loop_amp_orders']: 2376 for amp_number in SO[1]: 2377 ampSO_list[amp_number-1]=overall_so_basis.index(SO[0])+1 2378 2379 replace_dict['loopAmpSO'] = '\n'.join(self.format_integer_list( 2380 ampSO_list,'LOOPAMPORDERS')) 2381 ampSO_list=[-1]*sum(len(el[1]) for el in amps_orders['born_amp_orders']) 2382 for SO in amps_orders['born_amp_orders']: 2383 for amp_number in SO[1]: 2384 ampSO_list[amp_number-1]=overall_so_basis.index(SO[0])+1 2385 replace_dict['BornAmpSO'] = '\n'.join(self.format_integer_list( 2386 ampSO_list,'BORNAMPORDERS')) 2387 2388 # We then go to the TIR setup 2389 # The first entry is the CutTools, we make sure it is available 2390 looplibs_av=['.TRUE.'] 2391 # one should be careful about the order in the following 2392 for tir_lib in ['pjfry','iregi','golem']: 2393 looplibs_av.append('.TRUE.' if tir_lib in self.all_tir and \ 2394 self.tir_available_dict[tir_lib] else '.FALSE.') 2395 replace_dict['data_looplibs_av']=','.join(looplibs_av) 2396 2397 # Helicity offset convention 2398 # For a given helicity, the attached integer 'i' means 2399 # 'i' in ]-inf;-HELOFFSET[ -> Helicity is equal, up to a sign, 2400 # to helicity number abs(i+HELOFFSET) 2401 # 'i' == -HELOFFSET -> Helicity is analytically zero 2402 # 'i' in ]-HELOFFSET,inf[ -> Helicity is contributing with weight 'i'. 2403 # If it is zero, it is skipped. 2404 # Typically, the hel_offset is 10000 2405 replace_dict['hel_offset'] = 10000 2406 2407 # Extract overall denominator 2408 # Averaging initial state color, spin, and identical FS particles 2409 den_factor_line = self.get_den_factor_line(matrix_element) 2410 replace_dict['den_factor_line'] = den_factor_line 2411 2412 # When the user asks for the polarized matrix element we must 2413 # multiply back by the helicity averaging factor 2414 replace_dict['hel_avg_factor'] = matrix_element.get_hel_avg_factor() 2415 2416 if write_auxiliary_files: 2417 # Write out the color matrix 2418 (CMNum,CMDenom) = self.get_color_matrix(matrix_element) 2419 CMWriter=open(pjoin('..','MadLoop5_resources', 2420 '%(proc_prefix)sColorNumFactors.dat'%matrix_element.rep_dict),'w') 2421 for ColorLine in CMNum: 2422 CMWriter.write(' '.join(['%d'%C for C in ColorLine])+'\n') 2423 CMWriter.close() 2424 CMWriter=open(pjoin('..','MadLoop5_resources', 2425 '%(proc_prefix)sColorDenomFactors.dat'%matrix_element.rep_dict),'w') 2426 for ColorLine in CMDenom: 2427 CMWriter.write(' '.join(['%d'%C for C in ColorLine])+'\n') 2428 CMWriter.close() 2429 2430 # Write out the helicity configurations 2431 HelConfigs=matrix_element.get_helicity_matrix() 2432 HelConfigWriter=open(pjoin('..','MadLoop5_resources', 2433 '%(proc_prefix)sHelConfigs.dat'%matrix_element.rep_dict),'w') 2434 for HelConfig in HelConfigs: 2435 HelConfigWriter.write(' '.join(['%d'%H for H in HelConfig])+'\n') 2436 HelConfigWriter.close() 2437 2438 # Extract helas calls 2439 born_ct_helas_calls, uvct_helas_calls = \ 2440 fortran_model.get_born_ct_helas_calls(matrix_element, 2441 squared_orders=squared_orders,split_orders=split_orders) 2442 coef_construction, coef_merging = fortran_model.get_coef_construction_calls(\ 2443 matrix_element,group_loops=self.group_loops, 2444 squared_orders=squared_orders,split_orders=split_orders) 2445 2446 loop_CT_calls = fortran_model.get_loop_CT_calls(matrix_element,\ 2447 group_loops=self.group_loops, 2448 squared_orders=squared_orders, split_orders=split_orders) 2449 # The proc_prefix must be replaced 2450 coef_construction = [c % matrix_element.rep_dict for c 2451 in coef_construction] 2452 loop_CT_calls = [lc % matrix_element.rep_dict for lc in loop_CT_calls] 2453 2454 file = open(os.path.join(self.template_dir,\ 2455 'loop_matrix_standalone.inc')).read() 2456 2457 # Setup the contextual environment which is used in the splitting 2458 # functions below 2459 context = self.get_context(matrix_element) 2460 # Decide here wether we need to split the loop_matrix.f file or not. 2461 # 200 is reasonable but feel free to change it. 2462 if not noSplit and (len(matrix_element.get_all_amplitudes())>200): 2463 file=self.split_HELASCALLS(writer,replace_dict,\ 2464 'helas_calls_split.inc',file,born_ct_helas_calls,\ 2465 'born_ct_helas_calls','helas_calls_ampb', 2466 required_so_broadcaster = 'CT_REQ_SO_DONE', 2467 continue_label = 2000, context = context) 2468 file=self.split_HELASCALLS(writer,replace_dict,\ 2469 'helas_calls_split.inc',file,uvct_helas_calls,\ 2470 'uvct_helas_calls','helas_calls_uvct', 2471 required_so_broadcaster = 'UVCT_REQ_SO_DONE', 2472 continue_label = 3000, context=context) 2473 file=self.split_HELASCALLS(writer,replace_dict,\ 2474 'helas_calls_split.inc',file,coef_construction,\ 2475 'coef_construction','coef_construction', 2476 required_so_broadcaster = 'LOOP_REQ_SO_DONE', 2477 continue_label = 4000, context=context) 2478 else: 2479 replace_dict['born_ct_helas_calls']='\n'.join(born_ct_helas_calls) 2480 replace_dict['uvct_helas_calls']='\n'.join(uvct_helas_calls) 2481 replace_dict['coef_construction']='\n'.join(coef_construction) 2482 2483 # For loop induced processes, always split the loop_CT_calls because 2484 # they are used both in loop_matrix.f and mp_compute_loop_coef.f so 2485 # that it is quite nice to have them placed in the same subroutine. 2486 if not noSplit and ( (len(matrix_element.get_all_amplitudes())>200) or 2487 not matrix_element.get('processes')[0].get('has_born')): 2488 file=self.split_HELASCALLS(writer,replace_dict,\ 2489 'helas_calls_split.inc',file,loop_CT_calls,\ 2490 'loop_CT_calls','loop_CT_calls', 2491 required_so_broadcaster = 'CTCALL_REQ_SO_DONE', 2492 continue_label = 5000, context=context) 2493 else: 2494 replace_dict['loop_CT_calls']='\n'.join(loop_CT_calls) 2495 2496 # For loop induced processes, add the 'loop_CT_calls' entry to the 2497 # general_replace_dict so that it can be used by 2498 # write_mp_compute_loop_coefs later 2499 matrix_element.rep_dict['loop_CT_calls']=replace_dict['loop_CT_calls'] 2500 2501 replace_dict['coef_merging']='\n'.join(coef_merging) 2502 2503 file = file % replace_dict 2504 number_of_calls = len(filter(lambda call: call.find('CALL LOOP') != 0, \ 2505 loop_CT_calls)) 2506 if writer: 2507 # Write the file 2508 writer.writelines(file,context=context) 2509 return number_of_calls 2510 else: 2511 # Return it to be written along with the others 2512 return number_of_calls, file
2513 2514 #=============================================================================== 2515 # LoopProcessExporterFortranSA 2516 #===============================================================================
2517 -class LoopProcessExporterFortranMatchBox(LoopProcessOptimizedExporterFortranSA, 2518 export_v4.ProcessExporterFortranMatchBox):
2519 """Class to take care of exporting a set of loop matrix elements in the 2520 Fortran format.""" 2521 2522 default_opt = {'clean': False, 'complex_mass':False, 2523 'export_format':'madloop_matchbox', 'mp':True, 2524 'loop_dir':'', 'cuttools_dir':'', 2525 'fortran_compiler':'gfortran', 2526 'output_dependencies':'external', 2527 'sa_symmetry':True} 2528 2529 2530
2531 - def get_color_string_lines(self, matrix_element):
2532 """Return the color matrix definition lines for this matrix element. Split 2533 rows in chunks of size n.""" 2534 2535 return export_v4.ProcessExporterFortranMatchBox.get_color_string_lines(matrix_element)
2536 2537
2538 - def get_JAMP_lines(self, *args, **opts):
2539 """Adding leading color part of the colorflow""" 2540 2541 return export_v4.ProcessExporterFortranMatchBox.get_JAMP_lines(self, *args, **opts)
2542
2543 - def get_ME_identifier(self, matrix_element, group_number = None, group_elem_number = None):
2544 """ To not mix notations between borns and virtuals we call it here also MG5 """ 2545 return 'MG5_%d_'%matrix_element.get('processes')[0].get('id')
2546 2547 2548 #=============================================================================== 2549 # LoopInducedExporter 2550 #===============================================================================
2551 -class LoopInducedExporterME(LoopProcessOptimizedExporterFortranSA):
2552 """ A class to specify all the functions common to LoopInducedExporterMEGroup 2553 and LoopInducedExporterMENoGroup (but not relevant for the original 2554 Madevent exporters)""" 2555 2556 madloop_makefile_name = 'makefile_MadLoop' 2557 2558
2559 - def __init__(self, *args, **opts):
2560 """ Initialize the process, setting the proc characteristics.""" 2561 super(LoopInducedExporterME, self).__init__(*args, **opts) 2562 self.proc_characteristic['loop_induced'] = True
2563
2564 - def get_context(self,*args,**opts):
2565 """ Make sure that the contextual variable MadEventOutput is set to 2566 True for this exporter""" 2567 2568 context = super(LoopInducedExporterME,self).get_context(*args,**opts) 2569 context['MadEventOutput'] = True 2570 return context
2571 2572
2573 - def get_source_libraries_list(self):
2574 """ Returns the list of libraries to be compiling when compiling the 2575 SOURCE directory. It is different for loop_induced processes and 2576 also depends on the value of the 'output_dependencies' option""" 2577 2578 libraries_list = super(LoopInducedExporterME,self).\ 2579 get_source_libraries_list() 2580 2581 if self.dependencies=='internal': 2582 libraries_list.append('$(LIBDIR)libcts.$(libext)') 2583 libraries_list.append('$(LIBDIR)libiregi.$(libext)') 2584 2585 return libraries_list
2586 2593
2594 - def copy_v4template(self, *args, **opts):
2595 """Pick the right mother functions 2596 """ 2597 # Call specifically the necessary building functions for the mixed 2598 # template setup for both MadEvent and MadLoop standalone 2599 LoopProcessExporterFortranSA.loop_additional_template_setup(self, 2600 copy_Source_makefile=False) 2601 2602 LoopProcessOptimizedExporterFortranSA.\ 2603 loop_optimized_additional_template_setup(self)
2604 2605 2606 #=========================================================================== 2607 # Create jpeg diagrams, html pages,proc_card_mg5.dat and madevent.tar.gz 2608 #===========================================================================
2609 - def finalize_v4_directory(self, matrix_elements, history = "", makejpg = False, 2610 online = False, compiler='g77'):
2611 """Function to finalize v4 directory, for inheritance. 2612 """ 2613 2614 self.proc_characteristic['loop_induced'] = True
2615 2616 # This can be uncommented if one desires to have the MadLoop 2617 # initialization performed at the end of the output phase. 2618 # Alternatively, one can simply execute the command 'initMadLoop' in 2619 # the madevent interactive interface after the output. 2620 # from madgraph.interface.madevent_interface import MadLoopInitializer 2621 # MadLoopInitializer.init_MadLoop(self.dir_path, 2622 # subproc_prefix=self.SubProc_prefix, MG_options=None) 2623
2624 - def write_tir_cache_size_include(self, writer):
2625 """Write the file 'tir_cache_size.inc' which sets the size of the TIR 2626 cache the the user wishes to employ and the default value for it. 2627 This can have an impact on MadLoop speed when using stability checks 2628 but also impacts in a non-negligible way MadLoop's memory footprint. 2629 It is therefore important that the user can chose its size.""" 2630 2631 # In this case of MadLoop+MadEvent output, we set it to 2 because we 2632 # gain further speed increase with a TIR cache of size 2 due to the 2633 # the fact that we call MadLoop once per helicity configuration in this 2634 # case. 2635 tir_cach_size = "parameter(TIR_CACHE_SIZE=2)" 2636 writer.writelines(tir_cach_size)
2637
2638 - def write_matrix_element_v4(self, writer, matrix_element, fortran_model, 2639 proc_id = None, config_map = [], subproc_number = None):
2640 """ Write it the wrapper to call the ML5 subroutine in the library.""" 2641 2642 # Generating the MadEvent wrapping ME's routines 2643 if not matrix_element.get('processes') or \ 2644 not matrix_element.get('diagrams'): 2645 return 0 2646 2647 if not isinstance(writer, writers.FortranWriter): 2648 raise writers.FortranWriter.FortranWriterError(\ 2649 "writer not FortranWriter") 2650 2651 replace_dict = copy.copy(matrix_element.rep_dict) 2652 2653 # Extract version number and date from VERSION file 2654 info_lines = self.get_mg5_info_lines() 2655 replace_dict['info_lines'] = info_lines 2656 2657 # Extract process info lines 2658 process_lines = self.get_process_info_lines(matrix_element) 2659 replace_dict['process_lines'] = process_lines 2660 2661 # Set proc_id 2662 # It can be set to None when write_matrix_element_v4 is called without 2663 # grouping. In this case the subroutine SMATRIX should take an empty 2664 # suffix. 2665 if proc_id is None: 2666 replace_dict['proc_id'] = '' 2667 else: 2668 replace_dict['proc_id'] = proc_id 2669 2670 #set the average over the number of initial helicities 2671 replace_dict['hel_avg_factor'] = matrix_element.get_hel_avg_factor() 2672 2673 # Extract helicity lines 2674 helicity_lines = self.get_helicity_lines(matrix_element) 2675 replace_dict['helicity_lines'] = helicity_lines 2676 2677 2678 # Extract ndiags 2679 ndiags = len(matrix_element.get('diagrams')) 2680 replace_dict['ndiags'] = ndiags 2681 2682 # Set define_iconfigs_lines 2683 replace_dict['define_iconfigs_lines'] = \ 2684 """INTEGER MAPCONFIG(0:LMAXCONFIGS), ICONFIG 2685 COMMON/TO_MCONFIGS/MAPCONFIG, ICONFIG""" 2686 2687 if proc_id: 2688 # Set lines for subprocess group version 2689 # Set define_iconfigs_lines 2690 replace_dict['define_iconfigs_lines'] += \ 2691 """\nINTEGER SUBDIAG(MAXSPROC),IB(2) 2692 COMMON/TO_SUB_DIAG/SUBDIAG,IB""" 2693 # Set set_amp2_line 2694 replace_dict['configID_in_matrix'] = "SUBDIAG(%s)"%proc_id 2695 else: 2696 # Standard running 2697 # Set set_amp2_line 2698 replace_dict['configID_in_matrix'] = "MAPCONFIG(ICONFIG)" 2699 2700 # If group_numer 2701 replace_dict['ml_prefix'] = \ 2702 self.get_ME_identifier(matrix_element, subproc_number, proc_id) 2703 2704 # Extract ncolor 2705 ncolor = max(1, len(matrix_element.get('color_basis'))) 2706 replace_dict['ncolor'] = ncolor 2707 2708 n_tot_diags = len(matrix_element.get_loop_diagrams()) 2709 replace_dict['n_tot_diags'] = n_tot_diags 2710 2711 file = open(pjoin(_file_path, \ 2712 'iolibs/template_files/%s' % self.matrix_file)).read() 2713 file = file % replace_dict 2714 2715 # Write the file 2716 writer.writelines(file) 2717 2718 return 0, ncolor
2719
2720 - def get_amp2_lines(self, *args, **opts):
2721 """Make sure the function is implemented in the daughters""" 2722 2723 raise NotImplemented, 'The function get_amp2_lines must be called in '+\ 2724 ' the daugthers of LoopInducedExporterME'
2725 2726 #=============================================================================== 2727 # LoopInducedExporterMEGroup 2728 #===============================================================================
2729 -class LoopInducedExporterMEGroup(LoopInducedExporterME, 2730 export_v4.ProcessExporterFortranMEGroup):
2731 """Class to take care of exporting a set of grouped loop induced matrix 2732 elements""" 2733 2734 matrix_file = "matrix_loop_induced_madevent_group.inc" 2735 2741
2742 - def write_source_makefile(self, *args, **opts):
2743 """Pick the correct write_source_makefile function from 2744 ProcessExporterFortranMEGroup""" 2745 2746 export_v4.ProcessExporterFortranMEGroup.write_source_makefile(self, 2747 *args, **opts)
2748
2749 - def copy_v4template(self, *args, **opts):
2750 """Pick the right mother functions 2751 """ 2752 # Call specifically the necessary building functions for the mixed 2753 # template setup for both MadEvent and MadLoop standalone 2754 2755 # Start witht the MadEvent one 2756 export_v4.ProcessExporterFortranMEGroup.copy_v4template(self,*args,**opts) 2757 2758 # Then the MadLoop-standalone related one 2759 LoopInducedExporterME.copy_v4template(self, *args, **opts)
2760
2761 - def finalize_v4_directory(self, *args, **opts):
2762 """Pick the right mother functions 2763 """ 2764 # Call specifically what finalize_v4_directory must be used, so that the 2765 # MRO doesn't interfere. 2766 2767 self.proc_characteristic['loop_induced'] = True 2768 2769 export_v4.ProcessExporterFortranMEGroup.finalize_v4_directory( 2770 self,*args,**opts) 2771 2772 # And the finilize_v4 from LoopInducedExporterME which essentially takes 2773 # care of MadLoop virtuals initialization 2774 LoopInducedExporterME.finalize_v4_directory(self,*args,**opts)
2775
2776 - def generate_subprocess_directory_v4(self, subproc_group, 2777 fortran_model,group_number):
2778 """Generate the Pn directory for a subprocess group in MadEvent, 2779 including the necessary matrix_N.f files, configs.inc and various 2780 other helper files""" 2781 2782 # Generate the MadLoop files 2783 calls = 0 2784 matrix_elements = subproc_group.get('matrix_elements') 2785 for ime, matrix_element in enumerate(matrix_elements): 2786 calls += self.generate_loop_subprocess(matrix_element,fortran_model, 2787 group_number = group_number, proc_id = str(ime+1), 2788 # group_number = str(subproc_group.get('number')), proc_id = str(ime+1), 2789 config_map = subproc_group.get('diagram_maps')[ime]) 2790 2791 # Then generate the MadEvent files 2792 export_v4.ProcessExporterFortranMEGroup.generate_subprocess_directory_v4( 2793 self, subproc_group,fortran_model,group_number) 2794 2795 return calls
2796
2797 - def get_amp2_lines(self, matrix_element, replace_dict, config_map):
2798 """Return the various replacement dictionary inputs necessary for the 2799 multichanneling amp2 definition for the loop-induced MadEvent output. 2800 """ 2801 2802 if not config_map: 2803 raise MadGraph5Error, 'A multi-channeling configuration map is '+\ 2804 ' necessary for the MadEvent Loop-induced output with grouping.' 2805 2806 nexternal, ninitial = matrix_element.get_nexternal_ninitial() 2807 2808 ret_lines = [] 2809 # In this case, we need to sum up all amplitudes that have 2810 # identical topologies, as given by the config_map (which 2811 # gives the topology/config for each of the diagrams 2812 diagrams = matrix_element.get('diagrams') 2813 2814 # Note that we need to use AMP2 number corresponding to the first 2815 # diagram number used for that AMP2. 2816 # The dictionary below maps the config ID to this corresponding first 2817 # diagram number 2818 config_index_map = {} 2819 # For each diagram number, the dictionary below gives the config_id it 2820 # belongs to or 0 if it doesn't belong to any. 2821 loop_amp_ID_to_config = {} 2822 2823 # Combine the diagrams with identical topologies 2824 config_to_diag_dict = {} 2825 for idiag, diag in enumerate(matrix_element.get('diagrams')): 2826 try: 2827 config_to_diag_dict[config_map[idiag]].append(idiag) 2828 except KeyError: 2829 config_to_diag_dict[config_map[idiag]] = [idiag] 2830 2831 for config in sorted(config_to_diag_dict.keys()): 2832 config_index_map[config] = (config_to_diag_dict[config][0] + 1) 2833 2834 # First add the UV and R2 counterterm amplitudes of each selected 2835 # diagram for the multichannel config 2836 CT_amp_numbers = [a.get('number') for a in \ 2837 sum([diagrams[idiag].get_ct_amplitudes() for \ 2838 idiag in config_to_diag_dict[config]], [])] 2839 2840 for CT_amp_number in CT_amp_numbers: 2841 loop_amp_ID_to_config[CT_amp_number] = config 2842 2843 # Now add here the loop amplitudes. 2844 loop_amp_numbers = [a.get('amplitudes')[0].get('number') 2845 for a in sum([diagrams[idiag].get_loop_amplitudes() for \ 2846 idiag in config_to_diag_dict[config]], [])] 2847 2848 for loop_amp_number in loop_amp_numbers: 2849 loop_amp_ID_to_config[loop_amp_number] = config 2850 2851 # Notice that the config_id's are not necessarily sequential here, so 2852 # the size of the config_index_map array has to be the maximum over all 2853 # config_ids. 2854 # config_index_map should never be empty unless there was no diagram, 2855 # so the expression below is ok. 2856 n_configs = max(config_index_map.keys()) 2857 replace_dict['nmultichannel_configs'] = n_configs 2858 2859 # We must fill the empty entries of the map with the dummy amplitude 2860 # number 0. 2861 conf_list = [(config_index_map[i] if i in config_index_map else 0) \ 2862 for i in range(1,n_configs+1)] 2863 # Now the placeholder 'nmultichannels' refers to the number of 2864 # multi-channels which are contributing, so we must filter out zeros. 2865 replace_dict['nmultichannels'] = len([_ for _ in conf_list if _!=0]) 2866 2867 # Now write the amp2 related inputs in the replacement dictionary 2868 res_list = [] 2869 chunk_size = 6 2870 for k in xrange(0, len(conf_list), chunk_size): 2871 res_list.append("DATA (config_index_map(i),i=%6r,%6r) /%s/" % \ 2872 (k + 1, min(k + chunk_size, len(conf_list)), 2873 ','.join(["%6r" % i for i in conf_list[k:k + chunk_size]]))) 2874 2875 replace_dict['config_index_map_definition'] = '\n'.join(res_list) 2876 2877 res_list = [] 2878 n_loop_amps = max(loop_amp_ID_to_config.keys()) 2879 amp_list = [loop_amp_ID_to_config[i] for i in \ 2880 sorted(loop_amp_ID_to_config.keys()) if i!=0] 2881 chunk_size = 6 2882 for k in xrange(0, len(amp_list), chunk_size): 2883 res_list.append("DATA (CONFIG_MAP(i),i=%6r,%6r) /%s/" % \ 2884 (k + 1, min(k + chunk_size, len(amp_list)), 2885 ','.join(["%6r" % i for i in amp_list[k:k + chunk_size]]))) 2886 2887 replace_dict['config_map_definition'] = '\n'.join(res_list) 2888 2889 return
2890 2891 #=============================================================================== 2892 # LoopInducedExporterMENoGroup 2893 #===============================================================================
2894 -class LoopInducedExporterMENoGroup(LoopInducedExporterME, 2895 export_v4.ProcessExporterFortranME):
2896 """Class to take care of exporting a set of individual loop induced matrix 2897 elements""" 2898 2899 matrix_file = "matrix_loop_induced_madevent.inc" 2900 2906
2907 - def write_source_makefile(self, *args, **opts):
2908 """Pick the correct write_source_makefile function from 2909 ProcessExporterFortran""" 2910 2911 super(export_v4.ProcessExporterFortranME,self).\ 2912 write_source_makefile(*args, **opts)
2913
2914 - def copy_v4template(self, *args, **opts):
2915 """Pick the right mother functions 2916 """ 2917 # Call specifically the necessary building functions for the mixed 2918 # template setup for both MadEvent and MadLoop standalone 2919 2920 # Start witht the MadEvent one 2921 export_v4.ProcessExporterFortranME.copy_v4template(self,*args,**opts) 2922 2923 # Then the MadLoop-standalone related one 2924 LoopInducedExporterME.copy_v4template(self, *args, **opts)
2925
2926 - def finalize_v4_directory(self, *args, **opts):
2927 """Pick the right mother functions 2928 """ 2929 2930 self.proc_characteristic['loop_induced'] = True 2931 # Call specifically what finalize_v4_directory must be used, so that the 2932 # MRO doesn't interfere. 2933 export_v4.ProcessExporterFortranME.finalize_v4_directory( 2934 self,*args,**opts) 2935 2936 # And the finilize_v4 from LoopInducedExporterME which essentially takes 2937 # care of MadLoop virtuals initialization 2938 LoopInducedExporterME.finalize_v4_directory(self,*args,**opts)
2939
2940 - def generate_subprocess_directory_v4(self, matrix_element, fortran_model, me_number):
2941 """Generate the Pn directory for a subprocess group in MadEvent, 2942 including the necessary matrix_N.f files, configs.inc and various 2943 other helper files""" 2944 2945 # Then generate the MadLoop files 2946 calls = self.generate_loop_subprocess(matrix_element,fortran_model, 2947 group_number = me_number) 2948 2949 2950 # First generate the MadEvent files 2951 calls += export_v4.ProcessExporterFortranME.generate_subprocess_directory_v4( 2952 self, matrix_element, fortran_model, me_number) 2953 return calls
2954
2955 - def get_amp2_lines(self, matrix_element, replace_dict, config_map):
2956 """Return the amp2(i) = sum(amp for diag(i))^2 lines""" 2957 2958 if config_map: 2959 raise MadGraph5Error, 'A configuration map should not be specified'+\ 2960 ' for the Loop induced exporter without grouping.' 2961 2962 nexternal, ninitial = matrix_element.get_nexternal_ninitial() 2963 # Get minimum legs in a vertex 2964 vert_list = [max(diag.get_vertex_leg_numbers()) for diag in \ 2965 matrix_element.get('diagrams') if diag.get_vertex_leg_numbers()!=[]] 2966 minvert = min(vert_list) if vert_list!=[] else 0 2967 2968 # Note that we need to use AMP2 number corresponding to the first 2969 # diagram number used for that AMP2. 2970 # The dictionary below maps the config ID to this corresponding first 2971 # diagram number 2972 config_index_map = {} 2973 # For each diagram number, the dictionary below gives the config_id it 2974 # belongs to or 0 if it doesn't belong to any. 2975 loop_amp_ID_to_config = {} 2976 2977 n_configs = 0 2978 for idiag, diag in enumerate(matrix_element.get('diagrams')): 2979 # Ignore any diagrams with 4-particle vertices. 2980 use_for_multichanneling = True 2981 if diag.get_vertex_leg_numbers()!=[] and max(diag.get_vertex_leg_numbers()) > minvert: 2982 use_for_multichanneling = False 2983 curr_config = 0 2984 else: 2985 n_configs += 1 2986 curr_config = n_configs 2987 2988 if not use_for_multichanneling: 2989 if 0 not in config_index_map: 2990 config_index_map[0] = idiag + 1 2991 else: 2992 config_index_map[curr_config] = idiag + 1 2993 2994 CT_amps = [ a.get('number') for a in diag.get_ct_amplitudes()] 2995 for CT_amp in CT_amps: 2996 loop_amp_ID_to_config[CT_amp] = curr_config 2997 2998 Loop_amps = [a.get('amplitudes')[0].get('number') 2999 for a in diag.get_loop_amplitudes()] 3000 for Loop_amp in Loop_amps: 3001 loop_amp_ID_to_config[Loop_amp] = curr_config 3002 3003 # Now write the amp2 related inputs in the replacement dictionary 3004 n_configs = len([k for k in config_index_map.keys() if k!=0]) 3005 replace_dict['nmultichannel_configs'] = n_configs 3006 # Now the placeholder 'nmultichannels' refers to the number of 3007 # multi-channels which are contributing which, in the non-grouped case 3008 # is always equal to the total number of multi-channels. 3009 replace_dict['nmultichannels'] = n_configs 3010 3011 res_list = [] 3012 conf_list = [config_index_map[i] for i in sorted(config_index_map.keys()) 3013 if i!=0] 3014 chunk_size = 6 3015 for k in xrange(0, len(conf_list), chunk_size): 3016 res_list.append("DATA (config_index_map(i),i=%6r,%6r) /%s/" % \ 3017 (k + 1, min(k + chunk_size, len(conf_list)), 3018 ','.join(["%6r" % i for i in conf_list[k:k + chunk_size]]))) 3019 3020 replace_dict['config_index_map_definition'] = '\n'.join(res_list) 3021 3022 res_list = [] 3023 n_loop_amps = max(loop_amp_ID_to_config.keys()) 3024 amp_list = [loop_amp_ID_to_config[i] for i in \ 3025 sorted(loop_amp_ID_to_config.keys()) if i!=0] 3026 chunk_size = 6 3027 for k in xrange(0, len(amp_list), chunk_size): 3028 res_list.append("DATA (CONFIG_MAP(i),i=%6r,%6r) /%s/" % \ 3029 (k + 1, min(k + chunk_size, len(amp_list)), 3030 ','.join(["%6r" % i for i in amp_list[k:k + chunk_size]]))) 3031 3032 replace_dict['config_map_definition'] = '\n'.join(res_list)
3033