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.
Functions/Subroutines
modhiggsj Module Reference

Functions/Subroutines

subroutine, public evalamp_hj (p, res)
 
subroutine me2_ggg_tree (j1, j2, j3, sprod, me2)
 
subroutine me2_qbqg_tree (j1, j2, j3, sprod, me2)
 
subroutine spinoru2 (n, p, za, zb, s)
 

Function/Subroutine Documentation

◆ evalamp_hj()

subroutine, public modhiggsj::evalamp_hj ( real(dp), dimension(4,4), intent(in)  p,
real(dp), dimension(-5:5,-5:5), intent(out)  res 
)

Definition at line 12 of file mod_HiggsJ.F90.

12  real(dp), intent(in) :: p(4,4)
13  complex(dp) :: heftcoupl
14  real(dp), intent(out) :: res(-5:5,-5:5)
15  real(dp) :: sprod(4,4)
16  complex(dp) :: za(4,4), zb(4,4)
17  real(dp) :: restmp!, restmpid
18  integer :: i,j
19 
20  heftcoupl = gs*alphas/(6.0_dp * pi * vev)
21  res=zero
22 
23  call spinoru2(3,(/-p(:,1),-p(:,2),p(:,4)/),za,zb,sprod)
24 
25  call me2_ggg_tree(1,2,3,sprod,res(0,0))
26  res(0,0) = res(0,0)*avegg
27 
28  do i=1,5
29  call me2_qbqg_tree(2,3,1,sprod,res(0,i))
30  res(0,i) = -res(0,i)*aveqg
31  !print *,res(0,i)
32  enddo
33 
34  do i=1,5
35  call me2_qbqg_tree(1,3,2,sprod,res(i,0))
36  res(i,0) = -res(i,0)*aveqg
37  !print *,res(i,0)
38  enddo
39 
40  do i=1,5
41  call me2_qbqg_tree(3,2,1,sprod,res(0,-i))
42  res(0,-i) = -res(0,-i)*aveqg
43  !print *,res(0,-i)
44  enddo
45 
46  do i=1,5
47  call me2_qbqg_tree(1,3,2,sprod,res(-i,0))
48  res(-i,0) = -res(-i,0)*aveqg
49  !if(res(-i,0).le.0d0) then
50  !print *,res(-i,0)
51  !endif
52  enddo
53 
54  do i=1,5
55  call me2_qbqg_tree(1,2,3,sprod,res(-i,i))
56  res(-i,i) = res(-i,i) * aveqq
57 ! res(i,-i) = res(i,-i) + restmp * aveqq
58 ! res(0, i) = res(0, i) + restmp * aveqg
59 ! res(0,-i) = res(0,-i) + restmp * aveqg
60 ! res(i, 0) = res(i, 0) + restmp * aveqg
61 ! res(-i,0) = res(-i,0) + restmp * aveqg
62  enddo
63 
64  do i=1,5
65  call me2_qbqg_tree(2,1,3,sprod,res(i,-i))
66  res(i,-i) = res(i,-i)*aveqq
67  enddo
68 
69  res = res * (heftcoupl**2)
70 !if(res(0,0).gt.10000d0)then
71  !do i = -5,5
72  !do j = -5,5
73  !print *, res(i,j),i,j
74  !enddo
75  !enddo
76 ! print *, dsqrt(scr(p(:,1),p(:,1))), dsqrt(scr(p(:,2),p(:,2))), dsqrt(scr(p(:,3),p(:,3))), dsqrt(scr(p(:,4),p(:,4)))
77 ! print*, p
78 !endif
79 

◆ me2_ggg_tree()

subroutine modhiggsj::me2_ggg_tree ( integer, intent(in)  j1,
integer, intent(in)  j2,
integer, intent(in)  j3,
real(dp), dimension(3,3), intent(in)  sprod,
real(dp), intent(out)  me2 
)
private

Definition at line 87 of file mod_HiggsJ.F90.

87  integer, intent(in) :: j1,j2,j3
88  real(dp), intent(in) :: sprod(3,3)
89  real(dp), intent(out) :: me2
90  real(dp) :: s12, s13, s23, qsq
91  real(dp), parameter :: col_gg = 4.0d0 * ca**2 * cf
92 
93  me2 = zero
94 
95  s12 = sprod(j1,j2)
96  s13 = sprod(j1,j3)
97  s23 = sprod(j2,j3)
98 
99  qsq = s12 + s13 + s23
100 
101  me2 = (s12**4 + s13**4 + s23**4 + qsq**4)/s12/s13/s23
102 
103  me2 = me2 * two * col_gg
104 
105  return
106 

◆ me2_qbqg_tree()

subroutine modhiggsj::me2_qbqg_tree ( integer, intent(in)  j1,
integer, intent(in)  j2,
integer, intent(in)  j3,
real(dp), dimension(3,3), intent(in)  sprod,
real(dp), intent(out)  me2 
)
private

Definition at line 112 of file mod_HiggsJ.F90.

112  integer, intent(in) :: j1, j2, j3
113  real(dp), intent(in) :: sprod(3,3)
114  real(dp), intent(out) :: me2
115  real(dp) :: s12, s13,s23
116  real(dp), parameter :: col_qg = two * xn * cf
117 
118  me2 = zero
119 
120  s12 = sprod(j1,j2)
121  s13 = sprod(j1,j3)
122  s23 = sprod(j2,j3)
123 
124  me2 = two * (s13**2 + s23**2)/s12 * col_qg
125 !if(me2.le.0d0)then
126 !print *, s12,s13,s23
127 !endif
128  return
129 

◆ spinoru2()

subroutine modhiggsj::spinoru2 ( integer, intent(in)  n,
real(dp), dimension(4,n), intent(in)  p,
complex(dp), dimension(n,n), intent(out)  za,
complex(dp), dimension(n,n), intent(out)  zb,
real(dp), dimension(n,n), intent(out)  s 
)
private

Definition at line 137 of file mod_HiggsJ.F90.

137  implicit none
138  integer, intent(in) :: n
139  real(dp), intent(in) :: p(4,n)
140  complex(dp), intent(out) :: za(n,n), zb(n,n)
141  real(dp), intent(out) :: s(n,n)
142  integer :: i,j
143  complex(dp) :: c23(n), f(n)
144  real(dp) :: rt(n)
145 
146  !---if one of the vectors happens to be zero this routine fails.
147  do j=1,n
148  za(j,j)=czero
149  zb(j,j)=za(j,j)
150 
151  !-----positive energy case
152  if (p(1,j) .gt. zero) then
153  rt(j)=sqrt(p(2,j)+p(1,j))
154  c23(j)=cmplx(p(4,j),-p(3,j),kind=dp)
155  f(j)=(one,zero)
156  else
157  !-----negative energy case
158  rt(j)=sqrt(-p(1,j)-p(2,j))
159  c23(j)=cmplx(-p(4,j),p(3,j),kind=dp)
160  f(j)=ci
161  endif
162  enddo
163 
164  do i=2,n
165 
166  do j=1,i-1
167  s(i,j)=two*scr(p(:,i),p(:,j))
168  za(i,j)=f(i)*f(j) &
169  *(c23(i)*cmplx(rt(j)/rt(i),kind=dp)-c23(j)*cmplx(rt(i)/rt(j),kind=dp))
170 
171  if (abs(s(i,j)).lt.1d-5) then
172  zb(i,j)=-(f(i)*f(j))**2*conjg(za(i,j))
173  else
174  zb(i,j)=-cmplx(s(i,j),kind=dp)/za(i,j)
175  endif
176 
177  za(j,i)=-za(i,j)
178  zb(j,i)=-zb(i,j)
179  s(j,i)=s(i,j)
180 
181  enddo
182 
183  enddo
184 
185  return
186 
hto_units::one
real *8, parameter one
Definition: CALLING_cpHTO.f:184
hto_units::ci
real *8, dimension(1:2) ci
Definition: CALLING_cpHTO.f:189
hto_units::zero
real *8, parameter zero
Definition: CALLING_cpHTO.f:185