subroutine ZZ_HZZ_gs(p,msq) implicit none ************************************************************************ * Author: J. M. Campbell * * July, 2002. * * * * Weak Boson Fusion by Z-Z exchange only * * This routine calculates the dipole subtraction terms * * for the process: * * q(-p1) + q(-p2) --> H(p34) + q(p5) + q(p6) * * | * * | * * | * * ---> Z(e-(p3)+e+(p4))+Z(mu-(p5)+mu+(p6)) * * * ************************************************************************ include 'constants.f' include 'ptilde.f' include 'qqgg.f' integer j,k,nd double precision p(mxpart,4),msq(maxd,-nf:nf,-nf:nf) double precision . msq19_7(-nf:nf,-nf:nf),msq29_8(-nf:nf,-nf:nf), . msq18_2(-nf:nf,-nf:nf),msq29_1(-nf:nf,-nf:nf), . msq17_2(-nf:nf,-nf:nf),msq28_1(-nf:nf,-nf:nf), . sub19_7(4),sub29_8(4),sub79_1(4),sub89_2(4), . sub18_2(4),sub29_1(4), . sub17_2(4),sub28_1(4), . dummy(-nf:nf,-nf:nf),dummyv(-nf:nf,-nf:nf),dsubv external ZZ_HZZ,donothing_gvec ndmax=6 do j=-nf,nf do k=-nf,nf do nd=1,ndmax msq(nd,j,k)=0d0 enddo enddo enddo c---- calculate the dipoles: initial-final and final-initial c---- note that we do not require the gg dipoles, so the v-type c---- entries are left as dummies call dips(1,p,1,9,7,sub19_7,dsubv,msq19_7,dummyv, . ZZ_HZZ,donothing_gvec) call dips(1,p,7,9,1,sub79_1,dsubv,dummy,dummyv, . ZZ_HZZ,donothing_gvec) call dips(2,p,2,9,8,sub29_8,dsubv,msq29_8,dummyv, . ZZ_HZZ,donothing_gvec) call dips(2,p,8,9,2,sub89_2,dsubv,dummy,dummyv, . ZZ_HZZ,donothing_gvec) call dips(3,p,1,7,2,sub17_2,dsubv,msq17_2,dummyv, . ZZ_HZZ,donothing_gvec) call dips(4,p,2,8,1,sub28_1,dsubv,msq28_1,dummyv, . ZZ_HZZ,donothing_gvec) call dips(5,p,1,8,2,sub18_2,dsubv,msq18_2,dummyv, . ZZ_HZZ,donothing_gvec) call dips(6,p,2,9,1,sub29_1,dsubv,msq29_1,dummyv, . ZZ_HZZ,donothing_gvec) do j=-nf,nf do k=-nf,nf if ((j .eq. 0) .and. (k .eq. 0)) then goto 20 elseif ((j .ne. 0) .and. (k .ne. 0)) then if ((j .gt. 0) .and. (k .lt. 0)) then c--- q-qb case msq(1,j,k)=2d0*cf*(sub19_7(qq)+sub79_1(qq))*msq19_7(j,k) msq(2,j,k)=2d0*cf*(sub29_8(qq)+sub89_2(qq))*msq29_8(j,k) elseif ((j .lt. 0) .and. (k .gt. 0)) then c--- qb-q case msq(1,j,k)=2d0*cf*(sub19_7(qq)+sub79_1(qq))*msq19_7(j,k) msq(2,j,k)=2d0*cf*(sub29_8(qq)+sub89_2(qq))*msq29_8(j,k) else c--- q-q and qb-qb cases msq(1,j,k)=2d0*cf*(sub19_7(qq)+sub79_1(qq))*msq19_7(j,k) msq(2,j,k)=2d0*cf*(sub29_8(qq)+sub89_2(qq))*msq29_8(j,k) endif C---qg elseif ((j .gt. 0) .and. (k .eq. 0)) then msq(6,j,k)=2d0*tr*sub29_1(qg)*( . +msq29_1(j,+1)+msq29_1(j,+2)+msq29_1(j,+3) . +msq29_1(j,+4)+msq29_1(j,+5)) msq(4,j,k)=2d0*tr*sub28_1(qg)*( . +msq28_1(j,-1)+msq28_1(j,-2)+msq28_1(j,-3) . +msq28_1(j,-4)+msq28_1(j,-5)) C---qbg elseif ((j .lt. 0) .and. (k .eq. 0)) then msq(6,j,k)=2d0*tr*sub29_1(qg)*( . +msq29_1(j,-5)+msq29_1(j,-4)+msq29_1(j,-3) . +msq29_1(j,-2)+msq29_1(j,-1)) msq(4,j,k)=2d0*tr*sub28_1(qg)*( . +msq28_1(j,+1)+msq28_1(j,+2)+msq28_1(j,+3) . +msq28_1(j,+4)+msq28_1(j,+5)) C---gq elseif ((j .eq. 0) .and. (k .gt. 0)) then msq(5,j,k)=2d0*tr*sub18_2(qg)*( . +msq18_2(+5,k)+msq18_2(+4,k)+msq18_2(+3,k) . +msq18_2(+2,k)+msq18_2(+1,k)) msq(3,j,k)=2d0*tr*sub17_2(qg)*( . +msq17_2(-1,k)+msq17_2(-2,k)+msq17_2(-3,k) . +msq17_2(-4,k)+msq17_2(-5,k)) C---gqb elseif ((j .eq. 0) .and. (k .lt. 0)) then msq(5,j,k)=2d0*tr*sub18_2(qg)*( . +msq18_2(-5,k)+msq18_2(-4,k)+msq18_2(-3,k) . +msq18_2(-2,k)+msq18_2(-1,k)) msq(3,j,k)=2d0*tr*sub17_2(qg)*( . +msq17_2(+1,k)+msq17_2(+2,k)+msq17_2(+3,k) . +msq17_2(+4,k)+msq17_2(+5,k)) endif 20 continue enddo enddo return end