JHUGen MELA  v2.4.1
Matrix element calculations as used in JHUGen. MELA is an important tool that was used for the Higgs boson discovery and for precise measurements of its structure and interactions. Please see the website https://spin.pha.jhu.edu/ and papers cited there for more details, and kindly cite those papers when using this code.
Public Member Functions | List of all members
collier_coefs::g_cll Interface Reference

Public Member Functions

subroutine g_main_cll (G, Guv, p10, p21, p32, p43, p54, p65, p60, p20, p31, p42, p53, p64, p50, p61, p30, p41, p52, p63, p40, p51, p62, m02, m12, m22, m32, m42, m52, m62, rmax, Gerr, id_in, Gerr2)
 
subroutine g_arrays_cll (G, Guv, MomInv, masses2, rmax, Gerr, Gerr2)
 
subroutine g_list_cll (G, Guv, p10, p21, p32, p43, p54, p65, p60, p20, p31, p42, p53, p64, p50, p61, p30, p41, p52, p63, p40, p51, p62, m02, m12, m22, m32, m42, m52, m62, rmax, Gerr, Gerr2)
 
subroutine g_arrays_list_cll (G, Guv, MomInv, masses2, rmax, Gerr, Gerr2)
 

Detailed Description

Definition at line 70 of file collier_coefs.F90.

Member Function/Subroutine Documentation

◆ g_arrays_cll()

subroutine collier_coefs::g_cll::g_arrays_cll ( double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  G,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  Guv,
double complex, dimension(21), intent(in)  MomInv,
double complex, dimension(0:6), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Gerr,
double precision, dimension(0:rmax), intent(out), optional  Gerr2 
)

Definition at line 3242 of file collier_coefs.F90.

3242 
3243  integer, intent(in) :: rmax
3244  double complex, intent(in) :: MomInv(21), masses2(0:6)
3245  double complex, intent(out) :: G(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
3246  double complex, intent(out) :: Guv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
3247  double precision, optional, intent(out) :: Gerr(0:rmax),Gerr2(0:rmax)
3248  double precision :: Gerraux(0:rmax),Gerr2aux(0:rmax)
3249 
3250  if (present(gerr)) then
3251  if (present(gerr2)) then
3252  call g_main_cll(g,guv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
3253  mominv(7),mominv(8),mominv(9),mominv(10),mominv(11),mominv(12), &
3254  mominv(13),mominv(14),mominv(15),mominv(16),mominv(17),mominv(18), &
3255  mominv(19),mominv(20),mominv(21),masses2(0),masses2(1), &
3256  masses2(2),masses2(3),masses2(4),masses2(5),masses2(6),rmax,gerr,gerr2=gerr2)
3257  else
3258  call g_main_cll(g,guv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
3259  mominv(7),mominv(8),mominv(9),mominv(10),mominv(11),mominv(12), &
3260  mominv(13),mominv(14),mominv(15),mominv(16),mominv(17),mominv(18), &
3261  mominv(19),mominv(20),mominv(21),masses2(0),masses2(1), &
3262  masses2(2),masses2(3),masses2(4),masses2(5),masses2(6),rmax,gerr)
3263  end if
3264  else
3265  if (present(gerr2)) then
3266  call g_main_cll(g,guv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
3267  mominv(7),mominv(8),mominv(9),mominv(10),mominv(11),mominv(12), &
3268  mominv(13),mominv(14),mominv(15),mominv(16),mominv(17),mominv(18), &
3269  mominv(19),mominv(20),mominv(21),masses2(0),masses2(1), &
3270  masses2(2),masses2(3),masses2(4),masses2(5),masses2(6),rmax,gerraux,gerr2=gerr2)
3271  else
3272  call g_main_cll(g,guv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
3273  mominv(7),mominv(8),mominv(9),mominv(10),mominv(11),mominv(12), &
3274  mominv(13),mominv(14),mominv(15),mominv(16),mominv(17),mominv(18), &
3275  mominv(19),mominv(20),mominv(21),masses2(0),masses2(1), &
3276  masses2(2),masses2(3),masses2(4),masses2(5),masses2(6),rmax,gerraux)
3277  end if
3278  end if
3279 

◆ g_arrays_list_cll()

subroutine collier_coefs::g_cll::g_arrays_list_cll ( double complex, dimension(:), intent(out)  G,
double complex, dimension(:), intent(out)  Guv,
double complex, dimension(21), intent(in)  MomInv,
double complex, dimension(0:6), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Gerr,
double precision, dimension(0:rmax), intent(out), optional  Gerr2 
)

Definition at line 3402 of file collier_coefs.F90.

3402 
3403  integer, intent(in) :: rmax
3404  double complex, intent(in) :: MomInv(21), masses2(0:6)
3405  double complex, intent(out) :: G(:),Guv(:)
3406  double precision, optional, intent(out) :: Gerr(0:rmax),Gerr2(0:rmax)
3407  logical :: eflag
3408 
3409  if (7.gt.nmax_cll) then
3410  call seterrflag_cll(-10)
3411  call errout_cll('G_cll','Nmax_cll smaller 7',eflag,.true.)
3412  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
3413  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 7'
3414  call propagateerrflag_cll
3415  return
3416  end if
3417  if (rmax.gt.rmax_cll) then
3418  call seterrflag_cll(-10)
3419  call errout_cll('G_cll','argument rmax larger than rmax_cll',eflag,.true.)
3420  write(nerrout_cll,*) 'rmax =',rmax
3421  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
3422  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
3423  call propagateerrflag_cll
3424  return
3425  end if
3426 
3427  call g_arrays_list_checked_cll(g,guv,mominv,masses2,rmax,gerr,gerr2)
3428 

◆ g_list_cll()

subroutine collier_coefs::g_cll::g_list_cll ( double complex, dimension(:), intent(out)  G,
double complex, dimension(:), intent(out)  Guv,
double complex, intent(in)  p10,
double complex, intent(in)  p21,
double complex, intent(in)  p32,
double complex, intent(in)  p43,
double complex, intent(in)  p54,
double complex, intent(in)  p65,
double complex, intent(in)  p60,
double complex, intent(in)  p20,
double complex, intent(in)  p31,
double complex, intent(in)  p42,
double complex, intent(in)  p53,
double complex, intent(in)  p64,
double complex, intent(in)  p50,
double complex, intent(in)  p61,
double complex, intent(in)  p30,
double complex, intent(in)  p41,
double complex, intent(in)  p52,
double complex, intent(in)  p63,
double complex, intent(in)  p40,
double complex, intent(in)  p51,
double complex, intent(in)  p62,
double complex, intent(in)  m02,
double complex, intent(in)  m12,
double complex, intent(in)  m22,
double complex, intent(in)  m32,
double complex, intent(in)  m42,
double complex, intent(in)  m52,
double complex, intent(in)  m62,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Gerr,
double precision, dimension(0:rmax), intent(out), optional  Gerr2 
)

Definition at line 3296 of file collier_coefs.F90.

3296 
3297  integer, intent(in) :: rmax
3298  double complex, intent(in) :: p10,p21,p32,p43,p54,p65,p60,p20,p31,p42,p53
3299  double complex, intent(in) :: p64,p50,p61,p30,p41,p52,p63,p40,p51,p62
3300  double complex, intent(in) :: m02,m12,m22,m32,m42,m52,m62
3301  double complex, intent(out) :: G(:),Guv(:)
3302  double precision, optional, intent(out) :: Gerr(0:rmax),Gerr2(0:rmax)
3303  logical :: eflag
3304 
3305  if (7.gt.nmax_cll) then
3306  call seterrflag_cll(-10)
3307  call errout_cll('G_cll','Nmax_cll smaller 7',eflag,.true.)
3308  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
3309  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 7'
3310  call propagateerrflag_cll
3311  return
3312  end if
3313  if (rmax.gt.rmax_cll) then
3314  call seterrflag_cll(-10)
3315  call errout_cll('G_cll','argument rmax larger than rmax_cll',eflag,.true.)
3316  write(nerrout_cll,*) 'rmax =',rmax
3317  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
3318  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
3319  call propagateerrflag_cll
3320  return
3321  end if
3322 
3323  call g_list_checked_cll(g,guv,p10,p21,p32,p43,p54,p65,p60,p20,p31,p42,p53, &
3324  p64,p50,p61,p30,p41,p52,p63,p40,p51,p62, &
3325  m02,m12,m22,m32,m42,m52,m62,rmax,gerr,gerr2)
3326 

◆ g_main_cll()

subroutine collier_coefs::g_cll::g_main_cll ( double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  G,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  Guv,
double complex, intent(in)  p10,
double complex, intent(in)  p21,
double complex, intent(in)  p32,
double complex, intent(in)  p43,
double complex, intent(in)  p54,
double complex, intent(in)  p65,
double complex, intent(in)  p60,
double complex, intent(in)  p20,
double complex, intent(in)  p31,
double complex, intent(in)  p42,
double complex, intent(in)  p53,
double complex, intent(in)  p64,
double complex, intent(in)  p50,
double complex, intent(in)  p61,
double complex, intent(in)  p30,
double complex, intent(in)  p41,
double complex, intent(in)  p52,
double complex, intent(in)  p63,
double complex, intent(in)  p40,
double complex, intent(in)  p51,
double complex, intent(in)  p62,
double complex, intent(in)  m02,
double complex, intent(in)  m12,
double complex, intent(in)  m22,
double complex, intent(in)  m32,
double complex, intent(in)  m42,
double complex, intent(in)  m52,
double complex, intent(in)  m62,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Gerr,
integer, intent(in), optional  id_in,
double precision, dimension(0:rmax), intent(out), optional  Gerr2 
)

Definition at line 3000 of file collier_coefs.F90.

3000 
3001  integer, intent(in) :: rmax
3002  double complex, intent(in) :: p10,p21,p32,p43,p54,p65,p60,p20,p31,p42,p53
3003  double complex, intent(in) :: p64,p50,p61,p30,p41,p52,p63,p40,p51,p62
3004  double complex, intent(in) :: m02,m12,m22,m32,m42,m52,m62
3005  double complex, intent(out) :: G(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
3006  double complex, intent(out) :: Guv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
3007  double precision, optional, intent(out) :: Gerr(0:rmax),Gerr2(0:rmax)
3008  double precision :: Gerraux(0:rmax),Gerr2aux(0:rmax)
3009  double precision :: Gacc(0:rmax), Gacc2(0:rmax),norm,norm_coli,norm_dd
3010  integer, optional, intent(in) :: id_in
3011  double complex :: args(28)
3012  double complex :: elimcminf2
3013  integer :: errflag,id
3014  logical :: mflag,eflag
3015  integer :: r,n1,n2,n3,n4,n5,n6
3016 
3017  if (7.gt.nmax_cll) then
3018  call seterrflag_cll(-10)
3019  call errout_cll('G_cll','Nmax_cll smaller 7',eflag,.true.)
3020  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
3021  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 7'
3022  call propagateerrflag_cll
3023  return
3024  end if
3025  if (rmax.gt.rmax_cll) then
3026  call seterrflag_cll(-10)
3027  call errout_cll('G_cll','argument rmax larger than rmax_cll',eflag,.true.)
3028  write(nerrout_cll,*) 'rmax =',rmax
3029  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
3030  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
3031  call propagateerrflag_cll
3032  return
3033  end if
3034 
3035  mflag=.true.
3036  if (present(id_in)) then
3037  mflag=.false.
3038  id = id_in
3039  else
3040  id = 0
3041  end if
3042 
3043  if (mflag) then
3044  ! set ID of master call
3045  args(1) = p10
3046  args(2) = p21
3047  args(3) = p32
3048  args(4) = p43
3049  args(5) = p54
3050  args(6) = p65
3051  args(7) = p60
3052  args(8) = p20
3053  args(9) = p31
3054  args(10) = p42
3055  args(11) = p53
3056  args(12) = p64
3057  args(13) = p50
3058  args(14) = p61
3059  args(15) = p30
3060  args(16) = p41
3061  args(17) = p52
3062  args(18) = p63
3063  args(19) = p40
3064  args(20) = p51
3065  args(21) = p62
3066  args(22) = m02
3067  args(23) = m12
3068  args(24) = m22
3069  args(25) = m32
3070  args(26) = m42
3071  args(27) = m52
3072  args(28) = m62
3073  call setmasterfname_cll('G_cll')
3074  call setmastern_cll(7)
3075  call setmasterr_cll(rmax)
3076  call setmasterargs_cll(28,args)
3077 
3078  call settencache_cll(never_tenred_cll)
3079  end if
3080 
3081 
3082  select case (mode_cll)
3083 
3084  case (1)
3085  ! calculate loop integral using
3086  ! COLI implementation by AD/LH
3087 
3088  call calcg(g,guv,p10,p21,p32,p43,p54,p65,p60,p20,p31,p42,p53, &
3089  p64,p50,p61,p30,p41,p52,p63,p40,p51,p62, &
3090  m02,m12,m22,m32,m42,m52,m62,rmax,id,gerraux,gerr2aux)
3091 
3092  norm = abs(g(0,0,0,0,0,0,0))
3093  do r=1,rmax
3094  do n1=0,r
3095  do n2=0,r-n1
3096  do n3=0,r-n1-n2
3097  do n4=0,r-n1-n2-n3
3098  do n5=0,r-n1-n2-n3-n4
3099  n6=r-n1-n2-n3-n4-n5
3100  norm = max(norm,abs(g(0,n1,n2,n3,n4,n5,n6)))
3101  end do
3102  end do
3103  end do
3104  end do
3105  end do
3106  end do
3107  if (norm.eq.0d0) then
3108  norm = max(abs(p10),abs(p21),abs(p32),abs(p43),abs(p54), &
3109  abs(p65),abs(p60),abs(p20),abs(p31),abs(p42), &
3110  abs(p53),abs(p64),abs(p50),abs(p61),abs(p30), &
3111  abs(p41),abs(p52),abs(p63),abs(p40),abs(p51), &
3112  abs(p62),abs(m02),abs(m12),abs(m22),abs(m32), &
3113  abs(m42),abs(m52),abs(m62))
3114  if(norm.ne.0d0) then
3115  norm=1d0/norm**5
3116  else
3117  norm=1d0/muir2_cll**4
3118  end if
3119  end if
3120  gacc = gerraux/norm
3121  gacc2 = gerr2aux/norm
3122 
3123  if (present(gerr)) gerr = gerraux
3124  if (present(gerr2)) gerr2 = gerr2aux
3125 
3126  if (mflag) call propagateaccflag_cll(gacc,rmax)
3127 
3128 
3129  case (2)
3130  call seterrflag_cll(-10)
3131  call errout_cll('G_cll','7-point functions not implemented in DD library',eflag)
3132  if(eflag) then
3133  write(nerrout_cll,*) 'G_cll: 7-point functions not implemented in DD library'
3134  write(nerrout_cll,*) 'G_cll: --> use COLI implementation (mode_cll=1)'
3135  end if
3136 
3137 
3138  case (3)
3139  ! calculate loop integral using
3140  ! COLI implementation by AD/LH
3141 
3142  call calcg(g,guv,p10,p21,p32,p43,p54,p65,p60,p20,p31,p42,p53, &
3143  p64,p50,p61,p30,p41,p52,p63,p40,p51,p62, &
3144  m02,m12,m22,m32,m42,m52,m62,rmax,id,gerraux,gerr2aux)
3145 
3146  norm = abs(g(0,0,0,0,0,0,0))
3147  do r=1,rmax
3148  do n1=0,r
3149  do n2=0,r-n1
3150  do n3=0,r-n1-n2
3151  do n4=0,r-n1-n2-n3
3152  do n5=0,r-n1-n2-n3-n4
3153  n6=r-n1-n2-n3-n4-n5
3154  norm = max(norm,abs(g(0,n1,n2,n3,n4,n5,n6)))
3155  end do
3156  end do
3157  end do
3158  end do
3159  end do
3160  end do
3161  if (norm.eq.0d0) then
3162  norm = max(abs(p10),abs(p21),abs(p32),abs(p43),abs(p54), &
3163  abs(p65),abs(p60),abs(p20),abs(p31),abs(p42), &
3164  abs(p53),abs(p64),abs(p50),abs(p61),abs(p30), &
3165  abs(p41),abs(p52),abs(p63),abs(p40),abs(p51), &
3166  abs(p62),abs(m02),abs(m12),abs(m22),abs(m32), &
3167  abs(m42),abs(m52),abs(m62))
3168  if(norm.ne.0d0) then
3169  norm=1d0/norm**5
3170  else
3171  norm=1d0/muir2_cll**4
3172  end if
3173  end if
3174  gacc = gerraux/norm
3175  gacc2 = gerr2aux/norm
3176 
3177  if (present(gerr)) gerr = gerraux
3178  if (present(gerr2)) gerr2 = gerr2aux
3179 
3180  if (mflag) call propagateaccflag_cll(gacc,rmax)
3181 
3182  call seterrflag_cll(-10)
3183  call errout_cll('G_cll','7-point functions not implemented in DD library',eflag)
3184  if(eflag) then
3185  write(nerrout_cll,*) 'G_cll: 7-point functions not implemented in DD library'
3186  write(nerrout_cll,*) 'G_cll: --> use COLI implementation (mode_cll=1)'
3187  end if
3188 
3189 
3190  end select
3191 
3192  if (mflag) call propagateerrflag_cll
3193 
3194  if (monitoring) then
3195  pointscntg_cll = pointscntg_cll + 1
3196 
3197  if(maxval(gacc).gt.reqacc_cll) accpointscntg_cll = accpointscntg_cll + 1
3198 
3199  if(maxval(gacc).gt.critacc_cll) then
3200  critpointscntg_cll = critpointscntg_cll + 1
3201  if ( critpointscntg_cll.le.noutcritpointsmax_cll(7) ) then
3202  call critpointsout_cll('G_cll',0,maxval(gacc), critpointscntg_cll)
3203  if( critpointscntg_cll.eq.noutcritpointsmax_cll(7)) then
3204  write(ncheckout_cll,*) ' Further output of Critical Points for G_cll suppressed '
3205  write(nerrout_cll,*)
3206  endif
3207  end if
3208  end if
3209 
3210 #ifdef CritPoints2
3211  if (mode_cll.ne.2) then
3212  if(maxval(gacc2).gt.reqacc_cll) accpointscntg2_cll = accpointscntg2_cll + 1
3213 
3214  if(maxval(gacc2).gt.critacc_cll) then
3215  critpointscntg2_cll = critpointscntg2_cll + 1
3216  if ( critpointscntg2_cll.le.noutcritpointsmax_cll(7) ) then
3217  call critpointsout2_cll('G_cll',0,maxval(gacc2), critpointscntg2_cll)
3218  if( critpointscntg2_cll.eq.noutcritpointsmax_cll(7)) then
3219  write(ncpout2_cll,*) ' Further output of Critical Points for G_cll suppressed '
3220  write(ncpout2_cll,*)
3221  endif
3222  end if
3223  end if
3224  end if
3225 #endif
3226 
3227  end if
3228 
3229 

The documentation for this interface was generated from the following file:
endif
O0 g endif() string(TOLOWER "$
Definition: CMakeLists.txt:143