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