You can not select more than 25 topics
Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.
115 lines
3.2 KiB
115 lines
3.2 KiB
!
|
|
! SPDX-License-Identifier: GPL-3.0-or-later
|
|
! Copyright (C) 2019-2022, respective authors of MCFM.
|
|
!
|
|
subroutine qqb_trigam(p,msq)
|
|
implicit none
|
|
include 'types.f'
|
|
c***********************************************************************
|
|
c Author: C. Williams *
|
|
c March, 2013. *
|
|
c Lowest order matrix element squared, averaged over initial colors *
|
|
c and spins *
|
|
c q(-p1)+qbar(-p2) --> gam(p3) + gam(p4) + gam(p5) *
|
|
c***********************************************************************
|
|
|
|
include 'constants.f'
|
|
include 'nf.f'
|
|
include 'mxpart.f'
|
|
include 'zprods_decl.f'
|
|
include 'ewcouple.f'
|
|
include 'ewcharge.f'
|
|
real(dp):: p(mxpart,4),msq(-nf:nf,-nf:nf)
|
|
integer:: j,h1,h2,h3,h4
|
|
complex(dp):: qqb(2,2,2,2)
|
|
real(dp):: qqbsum
|
|
real(dp):: fac,statfac
|
|
parameter(statfac=one/6._dp)
|
|
|
|
c--- initialize matrix elements
|
|
msq(:,:)=0._dp
|
|
|
|
call spinoru(5,p,za,zb)
|
|
|
|
c--- fill qqb helicity amplitudes
|
|
call amp_lo_3gam(1,2,3,4,5,za,zb,qqb)
|
|
c--- note that summed, squared qbq amplitudes are identical
|
|
c call amp_lo_3gam(2,1,3,4,5,za,zb,qbq)
|
|
|
|
qqbsum=0._dp
|
|
do h1=1,2
|
|
do h2=1,2
|
|
do h3=1,2
|
|
do h4=1,2
|
|
qqbsum=qqbsum+real(qqb(h1,h2,h3,h4)*conjg(qqb(h1,h2,h3,h4)))
|
|
enddo
|
|
enddo
|
|
enddo
|
|
enddo
|
|
|
|
c--- overall factor except for photon charge which is applied below
|
|
fac=esq**3*xn*8._dp*aveqq*statfac
|
|
|
|
do j=-nf,nf
|
|
if (j /= 0) then
|
|
msq(j,-j)=fac*Q(j)**6*qqbsum
|
|
endif
|
|
enddo
|
|
|
|
return
|
|
end
|
|
|
|
|
|
subroutine amp_lo_3gam(p1,p2,p3,p4,p5,za,zb,amp)
|
|
implicit none
|
|
include 'types.f'
|
|
include 'constants.f'
|
|
include 'mxpart.f'
|
|
include 'zprods_decl.f'
|
|
integer:: p1,p2,p3,p4,p5
|
|
complex(dp):: amp(2,2,2,2),trigam
|
|
|
|
c======= amplitudes that are zero
|
|
amp(2,1,1,1)=czip
|
|
amp(2,2,2,2)=czip
|
|
amp(1,2,2,2)=czip
|
|
amp(1,1,1,1)=czip
|
|
|
|
c======= 3,4,5 symmetry
|
|
amp(1,2,2,1)=trigam(p1,p2,p3,p4,p5,za,zb)
|
|
amp(1,2,1,2)=trigam(p1,p2,p3,p5,p4,za,zb)
|
|
amp(1,1,2,2)=trigam(p1,p2,p4,p5,p3,za,zb)
|
|
|
|
c======= line reversal
|
|
amp(2,2,2,1)=trigam(p2,p1,p3,p4,p5,za,zb)
|
|
amp(2,2,1,2)=trigam(p2,p1,p3,p5,p4,za,zb)
|
|
amp(2,1,2,2)=trigam(p2,p1,p4,p5,p3,za,zb)
|
|
|
|
c======= conjugation
|
|
amp(2,1,1,2)=-trigam(p1,p2,p3,p4,p5,zb,za)
|
|
amp(2,1,2,1)=-trigam(p1,p2,p3,p5,p4,zb,za)
|
|
amp(2,2,1,1)=-trigam(p1,p2,p4,p5,p3,zb,za)
|
|
c======= conjugation and line reversal
|
|
amp(1,1,1,2)=-trigam(p2,p1,p3,p4,p5,zb,za)
|
|
amp(1,1,2,1)=-trigam(p2,p1,p3,p5,p4,zb,za)
|
|
amp(1,2,1,1)=-trigam(p2,p1,p4,p5,p3,zb,za)
|
|
|
|
return
|
|
end
|
|
|
|
|
|
function trigam(p1,p2,p3,p4,p5,za,zb)
|
|
implicit none
|
|
include 'types.f'
|
|
complex(dp):: trigam
|
|
include 'mxpart.f'
|
|
include 'zprods_decl.f'
|
|
integer:: p1,p2,p3,p4,p5
|
|
c----- amplitude for q(-,-p1),qb(+,-p2),gam(3,+),gam(4,+),gam(5,-)
|
|
c----- all momenta outgoing
|
|
|
|
trigam=za(p2,p1)*za(p1,p5)**2/
|
|
& (za(p1,p3)*za(p1,p4)*za(p2,p3)*za(p2,p4))
|
|
|
|
return
|
|
end
|