JHUGen MELA  JHUGen v7.5.6, MELA v2.4.2
Matrix element calculations as used in JHUGen.
collier_coefs::c_cll Interface Reference

Public Member Functions

subroutine c_main_cll (C, Cuv, p10, p21, p20, m02, m12, m22, rmax, Cerr, id_in, Cerr2)
 
subroutine c_arrays_cll (C, Cuv, MomInv, masses2, rmax, Cerr, Cerr2)
 
subroutine c_list_cll (C, Cuv, p10, p21, p20, m02, m12, m22, rmax, Cerr, Cerr2)
 
subroutine c_arrays_list_cll (C, Cuv, MomInv, masses2, rmax, Cerr, Cerr2)
 

Detailed Description

Definition at line 46 of file collier_coefs.F90.

Member Function/Subroutine Documentation

◆ c_arrays_cll()

subroutine collier_coefs::c_cll::c_arrays_cll ( double complex, dimension(0:rmax/2,0:rmax,0:rmax), intent(out)  C,
double complex, dimension(0:rmax/2,0:rmax,0:rmax), intent(out)  Cuv,
double complex, dimension(3), intent(in)  MomInv,
double complex, dimension(0:2), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Cerr,
double precision, dimension(0:rmax), intent(out), optional  Cerr2 
)

Definition at line 1002 of file collier_coefs.F90.

1002 
1003  integer, intent(in) :: rmax
1004  double complex, intent(in) :: MomInv(3), masses2(0:2)
1005  double complex, intent(out) :: Cuv(0:rmax/2,0:rmax,0:rmax)
1006  double complex, intent(out) :: C(0:rmax/2,0:rmax,0:rmax)
1007  double precision, optional, intent(out) :: Cerr(0:rmax),Cerr2(0:rmax)
1008  double precision :: Cerraux(0:rmax),Cerr2aux(0:rmax)
1009 
1010  if (present(cerr)) then
1011  if (present(cerr2)) then
1012  call c_main_cll(c,cuv,mominv(1),mominv(2),mominv(3), &
1013  masses2(0),masses2(1),masses2(2),rmax,cerr,cerr2=cerr2)
1014  else
1015  call c_main_cll(c,cuv,mominv(1),mominv(2),mominv(3), &
1016  masses2(0),masses2(1),masses2(2),rmax,cerr,cerr2=cerr2aux)
1017  end if
1018  else
1019  if (present(cerr2)) then
1020  call c_main_cll(c,cuv,mominv(1),mominv(2),mominv(3), &
1021  masses2(0),masses2(1),masses2(2),rmax,cerraux,cerr2=cerr2)
1022  else
1023  call c_main_cll(c,cuv,mominv(1),mominv(2),mominv(3), &
1024  masses2(0),masses2(1),masses2(2),rmax,cerraux,cerr2=cerr2aux)
1025  end if
1026  end if
1027 

◆ c_arrays_list_cll()

subroutine collier_coefs::c_cll::c_arrays_list_cll ( double complex, dimension(:), intent(out)  C,
double complex, dimension(:), intent(out)  Cuv,
double complex, dimension(3), intent(in)  MomInv,
double complex, dimension(0:2), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Cerr,
double precision, dimension(0:rmax), intent(out), optional  Cerr2 
)

Definition at line 1126 of file collier_coefs.F90.

1126 
1127  integer, intent(in) :: rmax
1128  double complex, intent(in) :: MomInv(3), masses2(0:2)
1129  double complex, intent(out) :: Cuv(:),C(:)
1130  double precision, optional, intent(out) :: Cerr(0:rmax),Cerr2(0:rmax)
1131  logical :: eflag
1132 
1133  if (3.gt.nmax_cll) then
1134  call seterrflag_cll(-10)
1135  call errout_cll('C_cll','Nmax_cll smaller 3',eflag,.true.)
1136  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
1137  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 3'
1138  call propagateerrflag_cll
1139  return
1140  end if
1141  if (rmax.gt.rmax_cll) then
1142  call seterrflag_cll(-10)
1143  call errout_cll('C_cll','argument rmax larger than rmax_cll',eflag,.true.)
1144  write(nerrout_cll,*) 'rmax =',rmax
1145  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
1146  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
1147  call propagateerrflag_cll
1148  return
1149  end if
1150 
1151  call c_arrays_list_checked_cll(c,cuv,mominv,masses2,rmax,cerr,cerr2)
1152 

◆ c_list_cll()

subroutine collier_coefs::c_cll::c_list_cll ( double complex, dimension(:), intent(out)  C,
double complex, dimension(:), intent(out)  Cuv,
double complex, intent(in)  p10,
double complex, intent(in)  p21,
double complex, intent(in)  p20,
double complex, intent(in)  m02,
double complex, intent(in)  m12,
double complex, intent(in)  m22,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Cerr,
double precision, dimension(0:rmax), intent(out), optional  Cerr2 
)

Definition at line 1040 of file collier_coefs.F90.

1040 
1041  integer, intent(in) :: rmax
1042  double complex, intent(in) :: p10,p21,p20,m02,m12,m22
1043  double complex, intent(out) :: Cuv(:),C(:)
1044  double precision, optional, intent(out) :: Cerr(0:rmax),Cerr2(0:rmax)
1045  logical :: eflag
1046 
1047  if (3.gt.nmax_cll) then
1048  call seterrflag_cll(-10)
1049  call errout_cll('C_cll','Nmax_cll smaller 3',eflag,.true.)
1050  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
1051  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 3'
1052  call propagateerrflag_cll
1053  return
1054  end if
1055  if (rmax.gt.rmax_cll) then
1056  call seterrflag_cll(-10)
1057  call errout_cll('C_cll','argument rmax larger than rmax_cll',eflag,.true.)
1058  write(nerrout_cll,*) 'rmax =',rmax
1059  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
1060  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
1061  call propagateerrflag_cll
1062  return
1063  end if
1064 
1065  call c_list_checked_cll(c,cuv,p10,p21,p20,m02,m12,m22,rmax,cerr,cerr2)
1066 

◆ c_main_cll()

subroutine collier_coefs::c_cll::c_main_cll ( double complex, dimension(0:rmax/2,0:rmax,0:rmax), intent(out)  C,
double complex, dimension(0:rmax/2,0:rmax,0:rmax), intent(out)  Cuv,
double complex, intent(in)  p10,
double complex, intent(in)  p21,
double complex, intent(in)  p20,
double complex, intent(in)  m02,
double complex, intent(in)  m12,
double complex, intent(in)  m22,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Cerr,
integer, intent(in), optional  id_in,
double precision, dimension(0:rmax), intent(out), optional  Cerr2 
)

Definition at line 699 of file collier_coefs.F90.

699 
700  integer, intent(in) :: rmax
701  double complex, intent(in) :: p10,p21,p20,m02,m12,m22
702  double precision :: q10,q21,q20
703  double complex :: mm02,mm12,mm22
704  double complex, intent(out) :: Cuv(0:rmax/2,0:rmax,0:rmax)
705  double complex, intent(out) :: C(0:rmax/2,0:rmax,0:rmax)
706  double precision, optional, intent(out) :: Cerr(0:rmax),Cerr2(0:rmax)
707  integer, optional, intent(in) :: id_in
708  double complex :: C2uv(0:rmax/2,0:rmax,0:rmax),C2(0:rmax/2,0:rmax,0:rmax)
709  double complex :: Ccoliuv(0:rmax,0:rmax,0:rmax),Ccoli(0:rmax,0:rmax,0:rmax)
710  double complex :: Cdduv(0:rmax,0:rmax,0:rmax)
711  double complex :: Cdd(0:rmax,0:rmax,0:rmax)
712  double precision :: Cerraux(0:rmax),Cerr2aux(0:rmax)
713  double complex :: elimcminf2
714  double complex args(6)
715  integer :: n0,rank,errflag,id
716  double precision :: accrelDD(0:rmax_DD),accabsDD(0:rmax_DD)
717  double precision :: accrel2DD(0:rmax_DD),accabs2DD(0:rmax_DD)
718  double precision :: Cacc(0:rmax),norm,norm_coli,norm_dd,Cacc2(0:rmax),Cdiff(0:rmax)
719  integer :: accflagDD,errflagDD,NDD,rankDD
720  logical :: mflag,eflag
721  integer :: r,n1,n2
722 
723  if (3.gt.nmax_cll) then
724  call seterrflag_cll(-10)
725  call errout_cll('C_cll','Nmax_cll smaller 3',eflag,.true.)
726  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
727  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 3'
728  call propagateerrflag_cll
729  return
730  end if
731  if (rmax.gt.rmax_cll) then
732  call seterrflag_cll(-10)
733  call errout_cll('C_cll','argument rmax larger than rmax_cll',eflag,.true.)
734  write(nerrout_cll,*) 'rmax =',rmax
735  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
736  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
737  call propagateerrflag_cll
738  return
739  end if
740 
741  mflag=.true.
742  if (present(id_in)) then
743  mflag=.false.
744  id = id_in
745  else
746  id = 0
747  end if
748 
749 
750  if (mflag) then
751  ! set ID of master call
752  args(1) = p10
753  args(2) = p21
754  args(3) = p20
755  args(4) = m02
756  args(5) = m12
757  args(6) = m22
758  call setmasterfname_cll('C_cll')
759  call setmastern_cll(3)
760  call setmasterr_cll(rmax)
761  call setmasterargs_cll(6,args)
762 
763  call settencache_cll(never_tenred_cll)
764  end if
765 
766 
767  select case (mode_cll)
768 
769  case (1)
770  ! calculate loop integral using
771  ! COLI implementation by AD/LH
772 
773  call calcc(ccoli,ccoliuv,p10,p21,p20,m02,m12,m22,rmax,id,cerraux,cerr2aux)
774 
775  norm = abs(ccoli(0,0,0))
776  do r=1,rmax
777  do n1=0,r
778  n2=r-n1
779  norm = max(norm,abs(ccoli(0,n1,n2)))
780  end do
781  end do
782  if (norm.eq.0d0) then
783  norm = max(abs(p10),abs(p21),abs(p20),abs(m02),abs(m12),abs(m22))
784  if(norm.ne.0d0) then
785  norm=1d0/norm
786  else
787  norm=1d0/muir2_cll
788  end if
789  end if
790 
791  if (norm.ne.0d0) then
792  cacc = cerraux/norm
793  cacc2 = cerr2aux/norm
794  else
795  cacc = 0d0
796  cacc2 = 0d0
797  end if
798 
799  if (present(cerr)) cerr = cerraux
800  if (present(cerr2)) cerr2 = cerr2aux
801 
802  if (mflag) call propagateaccflag_cll(cacc,rmax)
803 
804  c(0:rmax/2,0:rmax,0:rmax) = ccoli(0:rmax/2,0:rmax,0:rmax)
805  cuv(0:rmax/2,0:rmax,0:rmax) = ccoliuv(0:rmax/2,0:rmax,0:rmax)
806 
807 
808  case (2)
809  ! calculate loop integral using
810  ! DD implementation by SD
811 
812  id=0
813 
814  ! replace small masses by DD-identifiers
815  q10 = dreal(getminf2dd_cll(p10))
816  q21 = dreal(getminf2dd_cll(p21))
817  q20 = dreal(getminf2dd_cll(p20))
818  mm02 = getminf2dd_cll(m02)
819  mm12 = getminf2dd_cll(m12)
820  mm22 = getminf2dd_cll(m22)
821 
822  rank = rmax
823  call c_dd(cdd,cdduv,q10,q21,q20,mm02,mm12,mm22,rank,id)
824  c(0:rank/2,0:rank,0:rank) = cdd(0:rank/2,0:rank,0:rank)
825  cuv(0:rank/2,0:rank,0:rank) = cdduv(0:rank/2,0:rank,0:rank)
826 
827  call ddgetacc(accreldd,accabsdd,accrel2dd,accabs2dd,ndd,rankdd,accflagdd,errflagdd,id)
828  if(present(cerr)) cerr(0:rmax) = accabsdd(0:rmax)
829  if(present(cerr2)) cerr2(0:rmax) = accabs2dd(0:rmax)
830 
831  norm = abs(c(0,0,0))
832  do r=1,rmax
833  do n1=0,r
834  n2=r-n1
835  norm = max(norm,abs(c(0,n1,n2)))
836  end do
837  end do
838  if (norm.eq.0d0) then
839  norm = max(abs(p10),abs(p21),abs(p20),abs(m02),abs(m12),abs(m22))
840  if(norm.ne.0d0) then
841  norm=1d0/norm
842  else
843  norm=1d0/muir2_cll
844  end if
845  end if
846  if (norm.ne.0d0) then
847  cacc = accabsdd(0:rmax)/norm
848  cacc2 = accabs2dd(0:rmax)/norm
849  else
850  cacc(r) = 0d0
851  cacc2(r) = 0d0
852  end if
853  if (mflag) call propagateaccflag_cll(cacc,rmax)
854 
855 
856  case (3)
857  ! cross-check mode
858  ! compare results for loop integral
859  ! from COLI implementation by AD/LH and
860  ! from DD implementation by SD
861 
862  ! calculate loop integral using COLI
863  call calcc(ccoli,ccoliuv,p10,p21,p20,m02,m12,m22,rmax,id,cerraux,cerr2aux)
864 
865  c(0:rmax/2,0:rmax,0:rmax) = ccoli(0:rmax/2,0:rmax,0:rmax)
866  cuv(0:rmax/2,0:rmax,0:rmax) = ccoliuv(0:rmax/2,0:rmax,0:rmax)
867 
868 
869  ! calculate loop integral using DD
870 
871  id=0
872 
873  ! replace small masses by DD-identifiers
874  q10 = dreal(getminf2dd_cll(p10))
875  q21 = dreal(getminf2dd_cll(p21))
876  q20 = dreal(getminf2dd_cll(p20))
877  mm02 = getminf2dd_cll(m02)
878  mm12 = getminf2dd_cll(m12)
879  mm22 = getminf2dd_cll(m22)
880 
881  rank = rmax
882  call c_dd(cdd,cdduv,q10,q21,q20,mm02,mm12,mm22,rank,id)
883  c2(0:rank/2,0:rank,0:rank) = cdd(0:rank/2,0:rank,0:rank)
884  c2uv(0:rank/2,0:rank,0:rank) = cdduv(0:rank/2,0:rank,0:rank)
885 
886  call ddgetacc(accreldd,accabsdd,accrel2dd,accabs2dd,ndd,rankdd,accflagdd,errflagdd,id)
887 
888  ! cross-check
889  norm_coli = abs(c(0,0,0))
890  norm_dd = abs(c2(0,0,0))
891  do r=1,rmax
892  do n1=0,r
893  n2=r-n1
894  norm_coli = max(norm_coli,abs(c(0,n1,n2)))
895  norm_dd = max(norm_dd,abs(c(0,n1,n2)))
896  end do
897  end do
898  if (norm_coli.eq.0d0) then
899  norm_coli = max(abs(p10),abs(p21),abs(p20),abs(m02),abs(m12),abs(m22))
900  if(norm_coli.ne.0d0) then
901  norm_coli=1d0/norm_coli
902  else
903  norm_coli=1d0/muir2_cll
904  end if
905  end if
906  if (norm_dd.eq.0d0) then
907  norm_dd = max(abs(p10),abs(p21),abs(p20),abs(m02),abs(m12),abs(m22))
908  if(norm_dd.ne.0d0) then
909  norm_dd=1d0/norm_dd
910  else
911  norm_dd=1d0/muir2_cll
912  end if
913  end if
914  norm = min(norm_coli,norm_dd)
915 
916  call checkcoefsc_cll(c,c2,p10,p21,p20,m02,m12,m22,rmax,norm,cdiff)
917 
918 
919  if (cerraux(rmax).lt.accabsdd(rmax)) then
920  if (present(cerr)) cerr = max(cerraux,cdiff)
921  if (present(cerr2)) cerr2 = cerr2aux
922  cacc = max(cerraux/norm_coli,cdiff/norm)
923  cacc2 = cerr2aux/norm_coli
924  if (monitoring) pointscntc_coli = pointscntc_coli + 1
925  else
926  c = c2
927  cuv = c2uv
928  if (present(cerr)) cerr = max(accabsdd(0:rmax),cdiff)
929  if (present(cerr2)) cerr2 = accabs2dd(0:rmax)
930  cacc = max(accabsdd(0:rmax)/norm_dd,cdiff/norm)
931  cacc2 = accabs2dd(0:rmax)/norm_dd
932  if (monitoring) pointscntc_dd = pointscntc_dd + 1
933  end if
934 
935  if (mflag) call propagateaccflag_cll(cacc,rmax)
936 
937  end select
938 
939  if (mflag) call propagateerrflag_cll
940 
941  if (monitoring) then
942  pointscntc_cll = pointscntc_cll + 1
943 
944  if(maxval(cacc).gt.reqacc_cll) accpointscntc_cll = accpointscntc_cll + 1
945 
946  if(maxval(cacc).gt.critacc_cll) then
947  critpointscntc_cll = critpointscntc_cll + 1
948  if ( critpointscntc_cll.le.noutcritpointsmax_cll(3) ) then
949  call critpointsout_cll('C_cll',0,maxval(cacc), critpointscntc_cll)
950  if( critpointscntc_cll.eq.noutcritpointsmax_cll(3)) then
951  write(ncpout_cll,*) ' Further output of Critical Points for C_cll suppressed '
952  write(ncpout_cll,*)
953  endif
954  end if
955 
956 ! write(ncpout_cll,*) 'Cerr_coli =',Cerraux
957 ! write(ncpout_cll,*) 'Cerr_dd =',accabsDD(0:rmax)
958 ! write(ncpout_cll,*) 'norm_coli =',norm_coli
959 ! write(ncpout_cll,*) 'norm_dd =',norm_dd
960 ! write(ncpout_cll,*) 'Cacc_coli =',Cerraux/norm_coli
961 ! write(ncpout_cll,*) 'Cacc_dd =',accabsDD(0:rmax)/norm_dd
962 ! write(ncpout_cll,*) 'Cdiff =',Cdiff
963 ! write(ncpout_cll,*) 'norm =',norm
964 ! write(ncpout_cll,*) 'Caccdiff =',Cdiff/norm
965 !
966 ! write(ncpout_cll,*) 'Cacc =',Cacc
967 ! write(ncpout_cll,*) 'Cacc =',Cerraux(rmax).lt.accabsDD(rmax)
968 !
969 
970  end if
971 
972 #ifdef CritPoints2
973  if(maxval(cacc2).gt.reqacc_cll) accpointscntc2_cll = accpointscntc2_cll + 1
974 
975  if(maxval(cacc2).gt.critacc_cll) then
976  critpointscntc2_cll = critpointscntc2_cll + 1
977  if ( critpointscntc2_cll.le.noutcritpointsmax_cll(3) ) then
978  call critpointsout2_cll('C_cll',0,maxval(cacc2), critpointscntc2_cll)
979  if( critpointscntc2_cll.eq.noutcritpointsmax_cll(3)) then
980  write(ncpout2_cll,*) ' Further output of Critical Points for C_cll suppressed '
981  write(ncpout2_cll,*)
982  endif
983  end if
984  end if
985 #endif
986 
987  end if
988 
989 

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