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
hto_hbb_cp Module Reference

Functions/Subroutines

real *8 function hto_sshh (x)
 
real *8 function hto_shh (muhr, scal, rgh)
 

Function/Subroutine Documentation

◆ hto_shh()

real*8 function hto_hbb_cp::hto_shh ( real*8  muhr,
real*8  scal,
real*8  rgh 
)

Definition at line 3611 of file CALLING_cpHTO.f.

3611  USE hto_ferbos
3612  USE hto_a_cmplx
3613  USE hto_aux_hcp
3614  USE hto_aux_hbb
3615  USE hto_units
3616  USE hto_acmplx_pro
3617  USE hto_acmplx_rat
3618  USE hto_full_ln
3619  USE hto_ln_2_riemann
3620  USE hto_masses
3621  USE hto_set_phys_const
3622  USE hto_riemann
3623  USE hto_optcp
3624  USE hto_olas
3625  USE hto_common_niels
3626  USE hto_qcd
3627  USE hto_cmplx_rootz
3628 *
3629  IMPLICIT NONE
3630 *
3631  INTEGER nc,iz
3632  real*8 hto_shh,muhr,rgh,muhs,scal,scals,p2,xm0,str,sti,ewc,
3633  # sconv,asmur,emc,emb,emt,as_NLO,rmbs,rmcs,rmss,
3634  # crmbs,crmcs,lcxb,lcxc,lcxbs,lcxcs,lclxb,lclxc,qcdtop,neg
3635  real*8, dimension(2) :: axm0
3636  real*8, dimension(2,2) :: axms
3637  real*8, dimension(3,2,2) :: bmcb
3638  real*8, dimension(8,2,2) :: bmcf
3639  real*8, dimension(1,2,2) :: bmct
3640  real*8, dimension(3,2) :: bm0b,coefbb
3641  real*8, dimension(8,2) :: bm0f,coefbf
3642  real*8, dimension(1,2) :: bm0t,coefbt
3643  real*8, dimension(2) :: sh,shs,clh,b0sumb,b0sumf,cxp,ksumb,
3644  # ksumf,coefB1,coefB2,coefB3,coefB4,coefB5,coefB6,coefB7,
3645  # coefB8,coefB9,coefB10,coefB11,coefB12,ab2,ab3,
3646  # totalf,totalt,totalb,total,b0sumt,ksumt,xms,
3647  # b0part,cpt,ccxt,deltag,sww0,coefW1,coefW2,
3648  # ksumw,ccxts,cctsi,shi,clt,cxhw,cstsi,csmcts,cltmw,
3649  # b0sumw,sww,DW,b0sumw1,b0sumw2,cmxw,clmw,cxtw,
3650  # kfmsb,ktmsb,kbmsb,nloqcd,runbc,ttqcd,wccxt,wclt,wccxts,
3651  # wcxtw,wcltmw,ascal
3652 *
3653  INTERFACE
3654  SUBROUTINE hto_initalphas(asord,FR2,MUR,asmur,emc,emb,emt)
3655  USE hto_dzpar
3656  IMPLICIT NONE
3657  INTEGER asord
3658  real*8 fr2,mur,asmur,emc,emb,emt,hto_findalphasr0
3659  EXTERNAL hto_findalphasr0
3660  END SUBROUTINE hto_initalphas
3661  END INTERFACE
3662 *
3663  INTERFACE
3664  FUNCTION hto_alphas(MUR)
3665  USE hto_nffix
3666  USE hto_varflv
3667  USE hto_frrat
3668  USE hto_asinp
3669  USE hto_asfthr
3670  IMPLICIT NONE
3671  real*8 mur,hto_alphas
3672  END FUNCTION hto_alphas
3673  END INTERFACE
3674 *
3675  INTERFACE
3676  FUNCTION hto_quarkqcd(scal,psi,ps0i,xmsi,xm0i,type)
3677  # RESULT(value)
3678  USE hto_riemann
3679  USE hto_acmplx_pro
3680  USE hto_acmplx_rat
3681  USE hto_cmplx_root
3682  USE hto_cmplx_rootz
3683  USE hto_cmplx_srs_root
3684  USE hto_ln_2_riemann
3685  USE hto_full_ln
3686  USE hto_sp_fun
3687  USE hto_units
3688  IMPLICIT NONE
3689  INTEGER type
3690  real*8 scal,ps0i,xm0i
3691  real*8, dimension(2) :: value,psi,xmsi
3692  END FUNCTION hto_quarkqcd
3693  END INTERFACE
3694 *
3695  muhs= muhr*muhr
3696  scals= scal*scal
3697 *
3698  asmur= 0.12018d0
3699  emc= 1.4d0
3700  emb= 4.75d0
3701  emt= mt
3702  iz= 1
3703  CALL hto_initalphas(iz,one,mz,asmur,emc,emb,emt)
3704  als= hto_alphas(scal)
3705  as_nlo= als/pi
3706 *
3707  runbc= hto_run_bc(scal)
3708  crmbs= runbc(2)*runbc(2)
3709  crmcs= runbc(1)*runbc(1)
3710 *
3711  lcxb= crmbs/scals
3712  lcxc= crmcs/scals
3713 *
3714  lcxbs= lcxb*lcxb
3715  lcxcs= lcxc*lcxc
3716  lclxb= log(lcxb)
3717  lclxc= log(lcxc)
3718 *
3719  IF(gtop == 1) THEN
3720  imt= g_f/(8.d0*sqrt(2.d0)*pi)*(mt*mt-mw*mw)**2*
3721  # (mt*mt+2.d0*mw*mw)/mt**3
3722  ELSEIF(gtop == 0) THEN
3723  imt= 0.d0
3724  ELSE
3725  imt= yimt
3726  ENDIF
3727 *
3728  str= mt*mt-0.25d0*imt*imt
3729  sti= -mt*imt
3730  cpt(1)= str
3731  cpt(2)= sti
3732  wccxt= cpt/scals
3733  wccxts= wccxt.cp.wccxt
3734  ccxt= cpt/scals
3735  ccxts= ccxt.cp.ccxt
3736 *
3737  cctsi= co.cq.ccts
3738  cstsi= co.cq.csts
3739  csmcts= csts-ccts
3740  cmxw= -cxw
3741 *
3742  clmw= cmxw(1).fln.cmxw(2)
3743  clmw(2)= clmw(2)-2.d0*pi
3744 *
3745  sh(1)= muhs/scals
3746  sh(2)= -muhr*rgh/scals
3747 *
3748  cxhw= sh-cxw
3749  shs= sh.cp.sh
3750  shi= co.cq.sh
3751 *
3752  clh= sh(1).fln.sh(2)
3753  clt= ccxt(1).fln.ccxt(2)
3754  wclt= wccxt(1).fln.wccxt(2)
3755  cxtw= ccxt-cxw
3756  wcxtw= wccxt-cxw
3757  cltmw= cxtw(1).fln.cxtw(2)
3758  wcltmw= wcxtw(1).fln.wcxtw(2)
3759 *
3760 * w
3761  coefb1= (12.d0*cxw-4.d0*sh+(shs.cq.cxw))/64.d0
3762 *
3763 * z
3764  coefb2= (-4.d0*(sh.cq.ccts)+(shs.cq.cxw)+
3765  # 12.d0*(cxw.cq.cctq))/128.d0
3766 *
3767 * h
3768  coefb3= 9.d0/128.d0*(shs.cq.cxw)
3769 *
3770 * top
3771  coefb4= -3.d0/32.d0*((4.d0*ccxt-sh).cp.(ccxt.cq.cxw))
3772 *
3773 * light fermions
3774 *
3775  coefb5= -3.d0/32.d0*((4.d0*co*lcxb-sh).cq.cxw)*lcxb
3776 *
3777  coefb6= -1.d0/32.d0*((4.d0*co*cxtau-sh).cq.cxw)*cxtau
3778 *
3779  coefb7= -3.d0/32.d0*((4.d0*co*lcxc-sh).cq.cxw)*lcxc
3780 *
3781  coefb8= -1.d0/32.d0*((4.d0*co*cxmu-sh).cq.cxw)*cxmu
3782 *
3783  coefb9= -3.d0/32.d0*((4.d0*co*cxs-sh).cq.cxw)*cxs
3784 *
3785  coefb10= -3.d0/32.d0*((4.d0*co*cxd-sh).cq.cxw)*cxd
3786 *
3787  coefb11= -3.d0/32.d0*((4.d0*co*cxu-sh).cq.cxw)*cxu
3788 *
3789  coefb12= -1.d0/32.d0*((4.d0*co*cxe-sh).cq.cxw)*cxe
3790 *
3791  cxp(1)= muhs
3792  cxp(2)= -muhr*rgh
3793  p2= muhs
3794 *
3795  IF(ifb.eq.0) THEN
3796 *
3797  xms(1)= swr
3798  xms(2)= swi
3799  xm0= mw
3800  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3801  b0sumb= (coefb1.cp.b0part)
3802  xms(1)= szr
3803  xms(2)= szi
3804  xm0= mz
3805  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3806  b0sumb= b0sumb+(coefb2.cp.b0part)
3807  xms(1)= muhs
3808  xms(2)= -muhr*rgh
3809  xm0= muhr
3810  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3811  b0sumb= b0sumb+(coefb3.cp.b0part)
3812 *
3813  xms(1:2)= crmbs*co(1:2)
3814  xm0= sqrt(crmbs)
3815  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3816  b0sumf= (coefb5.cp.b0part)
3817  xms(1:2)= mtl*mtl*co(1:2)
3818  xm0= mtl
3819  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3820  b0sumf= b0sumf+(coefb6.cp.b0part)
3821  xms(1:2)= crmcs*co(1:2)
3822  xm0= sqrt(crmcs)
3823  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3824  b0sumf= b0sumf+(coefb7.cp.b0part)
3825  xms(1:2)= mm*mm*co(1:2)
3826  xm0= mm
3827  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3828  b0sumf= b0sumf+(coefb8.cp.b0part)
3829  xms(1:2)= msq*msq*co(1:2)
3830  xm0= msq
3831  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3832  b0sumf= b0sumf+(coefb9.cp.b0part)
3833  xms(1:2)= mdq*mdq*co(1:2)
3834  xm0= mdq
3835  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3836  b0sumf= b0sumf+(coefb10.cp.b0part)
3837  xms(1:2)= muq*muq*co(1:2)
3838  xm0= muq
3839  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3840  b0sumf= b0sumf+(coefb11.cp.b0part)
3841  xms(1:2)= me*me*co(1:2)
3842  xm0= me
3843  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3844  b0sumf= b0sumf+(coefb12.cp.b0part)
3845 *
3846  xms(1:2)= cpt(1:2)
3847  xm0= mt
3848  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3849  b0sumt= (coefb4.cp.b0part)
3850 *
3851  iz= 0
3852  xms(1:2)= crmbs*co(1:2)
3853  xm0= sqrt(crmbs)
3854  nloqcd= crmbs*hto_quarkqcd(scal,cxp,p2,xms,xm0,iz)
3855  xms(1:2)= crmcs*co(1:2)
3856  xm0= sqrt(crmcs)
3857  nloqcd= nloqcd+crmcs*hto_quarkqcd(scal,cxp,p2,xms,xm0,iz)
3858  IF(qcdc==0) nloqcd= 0.d0
3859 *
3860  xms(1:2)= cpt(1:2)
3861  xm0= mt
3862  iz= 1
3863  ttqcd= hto_quarkqcd(scal,cxp,p2,xms,xm0,iz)
3864  ttqcd= cpt.cp.ttqcd
3865  IF(qcdc==0) ttqcd= 0.d0
3866 *
3867  ELSEIF(ifb.eq.1) THEN
3868 *
3869  xms(1:2)= crmbs*co(1:2)
3870  xm0= sqrt(crmbs)
3871  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3872  b0sumf= (coefb5.cp.b0part)
3873  xms(1:2)= mtl*mtl*co(1:2)
3874  xm0= mtl
3875  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3876  b0sumf= b0sumf+(coefb6.cp.b0part)
3877  xms(1:2)= crmcs*co(1:2)
3878  xm0= sqrt(crmcs)
3879  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3880  b0sumf= b0sumf+(coefb7.cp.b0part)
3881  xms(1:2)= mm*mm*co(1:2)
3882  xm0= mm
3883  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3884  b0sumf= b0sumf+(coefb8.cp.b0part)
3885  xms(1:2)= msq*msq*co(1:2)
3886  xm0= msq
3887  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3888  b0sumf= b0sumf+(coefb9.cp.b0part)
3889  xms(1:2)= mdq*mdq*co(1:2)
3890  xm0= mdq
3891  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3892  b0sumf= b0sumf+(coefb10.cp.b0part)
3893  xms(1:2)= muq*muq*co(1:2)
3894  xm0= muq
3895  b0part=hto_b0af_em(scal,cxp,p2,xms,xm0)
3896  b0sumf= b0sumf+(coefb11.cp.b0part)
3897  xms(1:2)= me*me*co(1:2)
3898  xm0= me
3899  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3900  b0sumf= b0sumf+(coefb12.cp.b0part)
3901 *
3902  iz= 0
3903  xms(1:2)= crmbs*co(1:2)
3904  xm0= sqrt(crmbs)
3905  nloqcd= crmbs*hto_quarkqcd(scal,cxp,p2,xms,xm0,iz)
3906  xms(1:2)= crmcs*co(1:2)
3907  xm0= sqrt(crmcs)
3908  nloqcd= nloqcd+crmcs*hto_quarkqcd(scal,cxp,p2,xms,xm0,iz)
3909  IF(qcdc==0) nloqcd= 0.d0
3910 *
3911  ELSEIF(ifb.eq.2) THEN
3912 *
3913  xms(1:2)= cpt(1:2)
3914  xm0= mt
3915  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3916  b0sumt= (coefb4.cp.b0part)
3917 *
3918  xms(1:2)= cpt(1:2)
3919  xm0= mt
3920  iz= 1
3921  ttqcd= hto_quarkqcd(scal,cxp,p2,xms,xm0,iz)
3922  ttqcd= cpt.cp.ttqcd
3923  IF(qcdc==0) ttqcd= 0.d0
3924 *
3925  ELSEIF(ifb.eq.3) THEN
3926 *
3927  xms(1)= swr
3928  xms(2)= swi
3929  xm0= mw
3930  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3931  b0sumb= (coefb1.cp.b0part)
3932  xms(1)= szr
3933  xms(2)= szi
3934  xm0= mz
3935  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3936  b0sumb= b0sumb+(coefb2.cp.b0part)
3937  xms(1)= muhs
3938  xms(2)= -muhr*rgh
3939  xm0= muhr
3940  b0part= hto_b0af_em(scal,cxp,p2,xms,xm0)
3941  b0sumb= b0sumb+(coefb3.cp.b0part)
3942 *
3943  ENDIF
3944 *
3945  IF(ifb.eq.0.or.ifb.eq.1) THEN
3946  ksumf= -(sh.cq.cxw)*(cxe*clxe+cxmu*clxmu+
3947  # cxtau*clxtau+3.d0*(cxd*clxd+cxs*clxs+lcxb*lclxb+
3948  # cxu*clxu+lcxc*lclxc))/32.d0+
3949  # cxwi*(cxes+cxmus+cxtaus+3.d0*(
3950  # cxds+cxss+lcxbs+cxus+lcxcs))/8.d0
3951  ENDIF
3952  IF(ifb.eq.0.or.ifb.eq.2) THEN
3953  ksumt= 3.d0/32.d0*(4.d0*((ccxt.cq.cxw).cp.ccxt)
3954  # -(clt.cp.((sh.cq.cxw).cp.ccxt)))
3955  ENDIF
3956 *
3957  IF(ifb.eq.0.or.ifb.eq.3) THEN
3958 *
3959  ksumb= (-2.d0*((2.d0*co+cctsi+3.d0*(sh.cq.cxw)).cp.sh)
3960  # -(clcts.cp.((6.d0*cctsi-(sh.cq.cxw)).cp.sh))
3961  # +3.d0*(clw.cp.((4.d0*co+2.d0*cctsi-(sh.cq.cxw)).cp.sh))
3962  # -3.d0*(clh.cp.(shs.cq.cxw)))/128.d0
3963 *
3964  ENDIF
3965 *
3966  ewc= 4.d0*sqrt(2.d0)*g_f/pis
3967  sconv= muhs/scals
3968 *
3969  IF(ifb.eq.0.or.ifb.eq.1) THEN
3970  totalf= b0sumf+ksumf
3971  ENDIF
3972 *
3973  IF(ifb.eq.0.or.ifb.eq.2) THEN
3974  totalt= b0sumt+ksumt
3975 *
3976  neg= ewc*(swr*totalt(2)+swi*totalt(1))/sconv
3977  # -as_NLO*g_f/(sqrt(2.d0)*pis)*
3978  # (sh(1)*ttqcd(2)+sh(2)*ttqcd(1))/sconv
3979 *
3980  IF(neg < 0.d0) THEN
3981  totalt= 0.d0
3982  qcdtop= 0.d0
3983  pcnt= -1
3984  ELSE
3985  qcdtop= as_nlo*g_f/(sqrt(2.d0)*pis)*
3986  # (sh(1)*ttqcd(2)+sh(2)*ttqcd(1))/sconv
3987  pcnt= 1
3988  ENDIF
3989 *
3990  ENDIF
3991 *
3992  IF(ifb.eq.0.or.ifb.eq.3) THEN
3993  totalb= b0sumb+ksumb
3994  IF((swr*totalb(2)+swi*totalb(1)) < 0.d0) THEN
3995  totalb= 0.d0
3996  pcnt= -1
3997  ELSE
3998  pcnt= 1
3999  ENDIF
4000  ENDIF
4001 *
4002  total= totalf+totalt+totalb
4003 *
4004 *--- w self energies
4005 *
4006  deltag= 6.d0*co+0.5d0*(((7.d0*co-4.d0*csts).cq.csts).cp.clcts)
4007 *
4008  sww0= -(38.d0*cxw+6.d0*wccxt+7.d0*sh
4009  # -48.d0*(((wccxt.cq.sh).cq.cxw).cp.wccxt)+8.d0*(cxw.cq.sh))/128.d0
4010  # -3.d0/64.d0*((cxw-sh+(cxws.cq.cxhw)).cp.clh)
4011  # +3.d0/32.d0*(((co-4.d0*((wccxt.cq.sh).cq.cxw).cp.wccxt)).cp.wclt)
4012  # +((((8.d0*co-17.d0*cstsi+3.d0*cctsi).cp.cxw)
4013  # -6.d0*((cxw.cq.sh).cq.cctq)).cp.clcts)/64.d0
4014  # -((cxw.cq.sh).cq.cctq)/32.d0+5.d0/128.d0*(cxw.cq.ccts)
4015 *
4016  coefw1= -(((8.d0*co-(sh.cq.cxw)).cp.sh)*sh
4017  # -4.d0*((-12.d0*cxw+7.d0*sh).cp.cxw))/192.d0
4018 *
4019  coefw2= -((cxws.cq.csmcts).cp.(416.d0*co-192.d0*csts
4020  # -((132.d0*co-((12.d0*co+cctsi).cq.ccts)).cq.ccts)))/192.d0
4021 *
4022  cxp(1)= swr
4023  cxp(2)= swi
4024  p2= mw*mw
4025 *
4026  axms(1,1)= swr
4027  axms(1,2)= swi
4028  axm0(1)= mw
4029  axms(2,1)= muhs
4030  axms(2,2)= -muhr*rgh
4031  axm0(2)= muhr
4032  b0part= hto_b021_dm_cp(scal,cxp,p2,axms,axm0)
4033  b0sumw1= (coefw1.cp.b0part)
4034  b0sumw= (coefw1.cp.b0part)
4035 *
4036  axms(1,1)= szr
4037  axms(1,2)= szi
4038  axm0(1)= mz
4039  axms(2,1)= swr
4040  axms(2,2)= swi
4041  axm0(2)= mw
4042  b0part= hto_b021_dm_cp(scal,cxp,p2,axms,axm0)
4043  b0sumw2= (coefw2.cp.b0part)
4044  b0sumw= b0sumw+(coefw2.cp.b0part)
4045 *
4046  ksumw= -12.d0*((cxw
4047  # -0.5d0*((3.d0*co-(wccxts.cq.cxws)).cp.wccxt)).cp.wcltmw)
4048  # -((24.d0*cxw-((14.d0*co-(sh.cq.cxw)).cp.sh)).cp.clh)
4049  # +((36.d0*cxw-14.d0*sh-18.d0*((co-4.d0*(wccxt.cq.sh)).cp.wccxt)
4050  # +(shs.cq.cxw)).cp.clw)
4051  # -6.d0*(((2.d0*co+((cxwi-12.d0*shi).cp.wccxt)).cp.wccxt)
4052  # +1.d0/6.d0*((15.d0*co-(sh.cq.cxw)).cp.sh)
4053  # +2.d0/9.d0*((97.d0*co+9.d0*(cxw.cq.sh)).cp.cxw))
4054  # +(((cxw.cq.ccts).cp.(co-6.d0*(cxw.cq.sh))).cq.ccts)
4055  # -2.d0*(((cxw.cq.csmcts).cp.clcts).cp.(62.d0*co
4056  # -48.d0*csts-5.d0*cctsi))
4057  # -18.d0*(((cxws.cq.sh).cq.cctq).cp.clcts)
4058  # -72.d0*((wclt.cp.wccxts).cp.(shi-1.d0/12.d0*(wccxt.cq.cxws)))
4059  # +23.d0*(cxw.cq.ccts)
4060 *
4061  ksumw= ksumw/192.d0+3.d0/16.d0*(cxw.cp.(clw-clmw))
4062 *
4063  sww= b0sumw+ksumw
4064 *
4065  dw= (-sww+sww0)/sconv+deltag/16.d0
4066 * dw= 0.d0
4067 *
4068  IF(ifb.eq.0) THEN
4069  hto_shh= rgh/muhr*
4070  # (1.d0+EWC*(swr*DW(1)-swi*DW(2)))
4071  # -EWC*(swr*total(2)+swi*total(1))/sconv
4072  # +as_NLO*g_f/(sqrt(2.d0)*pis)*
4073  # (sh(1)*nloqcd(2)+sh(2)*nloqcd(1))/sconv
4074  # +qcdtop
4075 *
4076  ELSEIF(ifb.eq.1) THEN
4077  hto_shh= rgh/muhr*
4078  # (1.d0+EWC*(swr*DW(1)-swi*DW(2)))
4079  # -EWC*(swr*totalf(2)+swi*totalf(1))/sconv
4080  # +as_NLO*g_f/(sqrt(2.d0)*pis)*
4081  # (sh(1)*nloqcd(2)+sh(2)*nloqcd(1))/sconv
4082  ELSEIF(ifb.eq.2) THEN
4083  hto_shh= rgh/muhr*
4084  # (1.d0+EWC*(swr*DW(1)-swi*DW(2)))
4085  # -EWC*(swr*totalt(2)+swi*totalt(1))/sconv
4086  # +qcdtop
4087  ELSEIF(ifb.eq.3) THEN
4088  hto_shh= rgh/muhr*
4089  # (1.d0+EWC*(swr*DW(1)-swi*DW(2)))
4090  # -EWC*(swr*totalb(2)+swi*totalb(1))/sconv
4091  ENDIF
4092 *
4093  RETURN
4094 *

◆ hto_sshh()

real*8 function hto_hbb_cp::hto_sshh ( real*8, intent(in)  x)

Definition at line 3590 of file CALLING_cpHTO.f.

3590  USE hto_aux_hcp
3591 *
3592  IMPLICIT NONE
3593 *
3594 c I've commented the following line and added the other twos (Carlo Oleari)
3595 c REAL*8 HTO_SSHH,rgh,muc,scalc,x
3596  real*8 hto_sshh,rgh,muc,scalc
3597  real*8, INTENT(IN) :: x
3598 *
3599  muc= muhcp
3600  scalc= scalec
3601  rgh= x*muc
3602 *
3603  hto_sshh= hto_shh(muc,scalc,rgh)
3604 *
3605  RETURN
hto_masses::mm
real *8, parameter mm
Definition: CALLING_cpHTO.f:78
hto_acmplx_pro
Definition: CALLING_cpHTO.f:257
hto_masses::swr
real *8, parameter swr
Definition: CALLING_cpHTO.f:89
hto_aux_hbb::cxw
real *8, dimension(2) cxw
Definition: CALLING_cpHTO.f:227
hto_aux_hcp::pcnt
integer pcnt
Definition: CALLING_cpHTO.f:61
hto_aux_hcp::qcdc
integer qcdc
Definition: CALLING_cpHTO.f:61
hto_units::one
real *8, parameter one
Definition: CALLING_cpHTO.f:184
hto_masses::mz
real *8, parameter mz
Definition: CALLING_cpHTO.f:76
hto_cmplx_rootz
Definition: CALLING_cpHTO.f:680
hto_nffix
Definition: CALLING_cpHTO.f:2070
hto_masses::mdq
real *8, parameter mdq
Definition: CALLING_cpHTO.f:81
hto_set_phys_const
Definition: CALLING_cpHTO.f:239
hto_qcd
Definition: CALLING_cpHTO.f:3031
hto_aux_hcp::yimt
real *8 yimt
Definition: CALLING_cpHTO.f:62
hto_riemann
Definition: CALLING_cpHTO.f:195
hto_a_cmplx
Definition: CALLING_cpHTO.f:3356
hto_acmplx_rat::cq
real *8 function, dimension(2) cq(x, y)
Definition: CALLING_cpHTO.f:298
hto_units::co
real *8, dimension(1:2) co
Definition: CALLING_cpHTO.f:188
hto_set_phys_const::g_f
real *8 g_f
Definition: CALLING_cpHTO.f:240
hto_aux_hcp::cxmu
real *8 cxmu
Definition: CALLING_cpHTO.f:62
hto_acmplx_rat
Definition: CALLING_cpHTO.f:292
hto_masses::msq
real *8, parameter msq
Definition: CALLING_cpHTO.f:83
hto_aux_hcp::cxu
real *8 cxu
Definition: CALLING_cpHTO.f:62
hto_optcp
Definition: CALLING_cpHTO.f:233
hto_olas::hto_b0af_em
real *8 function, dimension(2) hto_b0af_em(scal, psi, ps0i, xmsi, xm0i)
Definition: CALLING_cpHTO.f:2800
hto_ln_2_riemann
Definition: CALLING_cpHTO.f:517
hto_asinp
Definition: CALLING_cpHTO.f:2061
hto_aux_hcp::cxs
real *8 cxs
Definition: CALLING_cpHTO.f:62
hto_masses::swi
real *8, parameter swi
Definition: CALLING_cpHTO.f:90
hto_masses::muq
real *8, parameter muq
Definition: CALLING_cpHTO.f:80
hto_aux_hcp
Definition: CALLING_cpHTO.f:60
hto_sp_fun
Definition: CALLING_cpHTO.f:753
hto_masses::mtl
real *8, parameter mtl
Definition: CALLING_cpHTO.f:79
hto_aux_hcp::clxe
real *8 clxe
Definition: CALLING_cpHTO.f:62
hto_findalphasr0
real *8 function hto_findalphasr0(ASI)
Definition: CALLING_cpHTO.f:2142
hto_units
Definition: CALLING_cpHTO.f:179
hto_dzpar
Definition: CALLING_cpHTO.f:2045
hto_full_ln
Definition: CALLING_cpHTO.f:464
hto_aux_hcp::cxe
real *8 cxe
Definition: CALLING_cpHTO.f:62
hto_a_cmplx::hto_b021_dm_cp
real *8 function, dimension(2) hto_b021_dm_cp(scal, psi, ps0i, xmsi, xm0i)
Definition: CALLING_cpHTO.f:3360
hto_aux_hcp::gtop
integer gtop
Definition: CALLING_cpHTO.f:61
hto_varflv
Definition: CALLING_cpHTO.f:2067
hto_ferbos
Definition: CALLING_cpHTO.f:219
hto_aux_hbb::cxws
real *8, dimension(2) cxws
Definition: CALLING_cpHTO.f:227
hto_full_ln::fln
real *8 function, dimension(2) fln(x, y)
Definition: CALLING_cpHTO.f:470
hto_aux_hcp::clcts
real *8, dimension(2) clcts
Definition: CALLING_cpHTO.f:67
hto_masses::mt
real *8 mt
Definition: CALLING_cpHTO.f:73
hto_aux_hcp::cxtau
real *8 cxtau
Definition: CALLING_cpHTO.f:62
hto_initalphas
subroutine hto_initalphas(IORD, FR2, MUR, ASMUR, MC, MB, MT)
Definition: CALLING_cpHTO.f:2089
hto_masses
Definition: CALLING_cpHTO.f:72
hto_masses::mw
real *8, parameter mw
Definition: CALLING_cpHTO.f:75
hto_ferbos::ifb
integer ifb
Definition: CALLING_cpHTO.f:220
hto_aux_hcp::scalec
real *8 scalec
Definition: CALLING_cpHTO.f:62
hto_olas
Definition: CALLING_cpHTO.f:2796
hto_aux_hbb
Definition: CALLING_cpHTO.f:225
hto_quarkqcd
real *8 function, dimension(2) hto_quarkqcd(scal, psi, ps0i, xmsi, xm0i, type)
Definition: CALLING_cpHTO.f:2851
hto_riemann::pis
real *8, parameter pis
Definition: CALLING_cpHTO.f:197
hto_aux_hcp::clxmu
real *8 clxmu
Definition: CALLING_cpHTO.f:62
hto_aux_hbb::clw
real *8, dimension(2) clw
Definition: CALLING_cpHTO.f:227
hto_asfthr
Definition: CALLING_cpHTO.f:2076
hto_masses::szi
real *8, parameter szi
Definition: CALLING_cpHTO.f:92
hto_cmplx_root
Definition: CALLING_cpHTO.f:659
hto_qcd::hto_run_bc
real *8 function, dimension(2) hto_run_bc(scal)
Definition: CALLING_cpHTO.f:3157
hto_aux_hcp::cxd
real *8 cxd
Definition: CALLING_cpHTO.f:62
hto_aux_hbb::ccts
real *8, dimension(2) ccts
Definition: CALLING_cpHTO.f:227
hto_aux_hcp::muhcp
real *8 muhcp
Definition: CALLING_cpHTO.f:62
hto_masses::szr
real *8, parameter szr
Definition: CALLING_cpHTO.f:91
hto_set_phys_const::imt
real *8 imt
Definition: CALLING_cpHTO.f:240
hto_aux_hbb::csts
real *8, dimension(2) csts
Definition: CALLING_cpHTO.f:227
hto_frrat
Definition: CALLING_cpHTO.f:2073
hto_set_phys_const::als
real *8 als
Definition: CALLING_cpHTO.f:240
hto_acmplx_pro::cp
real *8 function, dimension(2) cp(x, y)
Definition: CALLING_cpHTO.f:263
hto_masses::me
real *8, parameter me
Definition: CALLING_cpHTO.f:77
hto_common_niels
Definition: CALLING_cpHTO.f:717
hto_alphas
real *8 function hto_alphas(MUR)
Definition: CALLING_cpHTO.f:2269
hto_cmplx_srs_root
Definition: CALLING_cpHTO.f:618