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

Public Member Functions

subroutine bten_main_cll (TB, TBuv, MomVec, MomInv, masses2, rmax, TBerr)
 
subroutine bten_list_cll (TB, TBuv, MomVec, MomInv, masses2, rmax, TBerr)
 
subroutine bten_args_cll (TB, TBuv, p1vec, p10, m02, m12, rmax, TBerr)
 
subroutine bten_args_list_cll (TB, TBuv, p1vec, p10, m02, m12, rmax, TBerr)
 

Detailed Description

Definition at line 50 of file collier_tensors.F90.

Member Function/Subroutine Documentation

◆ bten_args_cll()

subroutine collier_tensors::bten_cll::bten_args_cll ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TB,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TBuv,
double complex, dimension(0:3), intent(in)  p1vec,
double complex, intent(in)  p10,
double complex, intent(in)  m02,
double complex, intent(in)  m12,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  TBerr 
)

Definition at line 1136 of file collier_tensors.F90.

1136 
1137  integer, intent(in) :: rmax
1138  double complex, intent(in) :: p1vec(0:3)
1139  double complex, intent(in) :: p10,m02,m12
1140  double complex, intent(out) :: TB(0:rmax,0:rmax,0:rmax,0:rmax)
1141  double complex, intent(out) :: TBuv(0:rmax,0:rmax,0:rmax,0:rmax)
1142  double precision, intent(out), optional :: TBerr(0:rmax)
1143  double complex :: TB2(0:rmax,0:rmax,0:rmax,0:rmax), TBuv2(0:rmax,0:rmax,0:rmax,0:rmax)
1144  double complex :: masses2(0:1),MomInv(1)
1145  double complex :: CB(0:rmax/2,0:rmax), CBuv(0:rmax/2,0:rmax)
1146  double precision :: CBerr(0:rmax),TBerr_aux(0:rmax),TBerr_aux2(0:rmax)
1147  double complex :: args(7)
1148  double precision :: TBdiff(0:rmax),norm(0:rmax),norm_coli,norm_dd,TBacc(0:rmax)
1149  integer :: r,n0,n1,n2,n3
1150  logical :: eflag
1151 
1152  if (2.gt.nmax_cll) then
1153  call seterrflag_cll(-10)
1154  call errout_cll('Bten_cll','Nmax_cll smaller 2',eflag,.true.)
1155  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
1156  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 2'
1157  call propagateerrflag_cll
1158  return
1159  end if
1160  if (rmax.gt.rmax_cll) then
1161  call seterrflag_cll(-10)
1162  call errout_cll('Bten_cll','argument rmax larger than rmax_cll',eflag,.true.)
1163  write(nerrout_cll,*) 'rmax =',rmax
1164  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
1165  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
1166  call propagateerrflag_cll
1167  return
1168  end if
1169 
1170  masses2(0) = m02
1171  masses2(1) = m12
1172  mominv(1) = p10
1173 
1174  ! set ID of master call
1175  args(1:4) = p1vec(0:)
1176  args(5) = p10
1177  args(6:7) = masses2(0:)
1178  call setmasterfname_cll('Bten_cll')
1179  call setmastern_cll(2)
1180  call setmasterr_cll(rmax)
1181  call setmasterargs_cll(7,args)
1182 
1183  call settencache_cll(tenred_cll-1)
1184 
1185  if (mode_cll.eq.3) then
1186  ! calculate tensor with coefficients from COLI
1187  mode_cll = 1
1188  call b_main_cll(cb,cbuv,p10,m02,m12,rmax,cberr,0)
1189  call calctensorb(tb,tbuv,tberr_aux,cb,cbuv,cberr,p1vec,rmax)
1190 
1191  ! calculate tensor with coefficients from DD
1192  mode_cll = 2
1193  call b_main_cll(cb,cbuv,p10,m02,m12,rmax,cberr,0)
1194  call calctensorb(tb2,tbuv2,tberr_aux2,cb,cbuv,cberr,p1vec,rmax)
1195 
1196  ! comparison --> take better result
1197  mode_cll = 3
1198  do r=0,rmax
1199  norm_coli=0d0
1200  norm_dd=0d0
1201  do n0=0,r
1202  do n1=0,r-n0
1203  do n2=0,r-n0-n1
1204  n3=r-n0-n1-n2
1205  norm_coli = max(norm_coli,abs(tb(n0,n1,n2,n3)))
1206  norm_dd = max(norm_dd,abs(tb2(n0,n1,n2,n3)))
1207  end do
1208  end do
1209  end do
1210  if (norm_coli.eq.0d0) then
1211  norm_coli = max(abs(mominv(1)),maxval(abs(masses2(0:1))))
1212  if(norm_coli.ne.0d0) then
1213  norm_coli=norm_coli**(real(r)/2)
1214  else
1215  norm_coli=muir2_cll**(real(r)/2)
1216  end if
1217  end if
1218  if (norm_dd.eq.0d0) then
1219  norm_dd = max(abs(mominv(1)),maxval(abs(masses2(0:1))))
1220  if(norm_dd.ne.0d0) then
1221  norm_dd=norm_dd**(real(r)/2)
1222  else
1223  norm_dd=muir2_cll**(real(r)/2)
1224  end if
1225  end if
1226  norm(r) = min(norm_coli,norm_dd)
1227  end do
1228 
1229  call checktensors_cll(tb,tb2,p1vec,mominv,masses2,norm,2,rmax,tbdiff)
1230 
1231  if (tberr_aux(rmax).lt.tberr_aux2(rmax)) then
1232  if (present(tberr)) tberr = max(tberr_aux,tbdiff*norm)
1233  do r=0,rmax
1234  tbacc(r) = max(tberr_aux(r)/norm(r),tbdiff(r))
1235  end do
1236  if (monitoring) pointscntbten_coli = pointscntbten_coli + 1
1237  else
1238  tb = tb2
1239  tbuv = tbuv2
1240  if (present(tberr)) tberr = max(tberr_aux2,tbdiff*norm)
1241  do r=0,rmax
1242  tbacc(r) = max(tberr_aux2(r)/norm(r),tbdiff(r))
1243  end do
1244  if (monitoring) pointscntbten_dd = pointscntbten_dd + 1
1245  end if
1246 
1247  else
1248  call b_main_cll(cb,cbuv,p10,m02,m12,rmax,cberr,0)
1249  call calctensorb(tb,tbuv,tberr_aux,cb,cbuv,cberr,p1vec,rmax)
1250  if (present(tberr)) tberr = tberr_aux
1251  norm = 0d0
1252  do r=0,rmax
1253  do n0=0,r
1254  do n1=0,r-n0
1255  do n2=0,r-n0-n1
1256  n3=r-n0-n1-n2
1257  norm(r) = max(norm(r),abs(tb(n0,n1,n2,n3)))
1258  end do
1259  end do
1260  end do
1261  if (norm(r).eq.0d0) then
1262  norm(r) = max(abs(mominv(1)),maxval(abs(masses2(0:1))))
1263  if(norm(r).ne.0d0) then
1264  norm(r)=norm(r)**(real(r)/2)
1265  else
1266  norm(r)=muir2_cll**(real(r)/2)
1267  end if
1268  end if
1269  tbacc(r) = tberr_aux(r)/norm(r)
1270  end do
1271 
1272  end if
1273 
1274  call propagateaccflag_cll(tbacc,rmax)
1275  call propagateerrflag_cll
1276 
1277  if (monitoring) then
1278  pointscntbten_cll = pointscntbten_cll + 1
1279 
1280  if(maxval(tbacc).gt.reqacc_cll) accpointscntbten_cll = accpointscntbten_cll + 1
1281 
1282  if(maxval(tbacc).gt.critacc_cll) then
1283  critpointscntbten_cll = critpointscntbten_cll + 1
1284  if ( critpointscntbten_cll.le.noutcritpointsmax_cll(2) ) then
1285  call critpointsout_cll('TBten_cll',0,maxval(tbacc),critpointscntbten_cll)
1286  if( critpointscntbten_cll.eq.noutcritpointsmax_cll(2)) then
1287  write(ncpout_cll,*) ' Further output of Critical Points for TBten_cll suppressed'
1288  write(ncpout_cll,*)
1289  endif
1290 #ifdef CritPoints2
1291  call critpointsout2_cll('TBten_cll',0,maxval(tbacc),critpointscntbten_cll)
1292  if( critpointscntbten_cll.eq.noutcritpointsmax_cll(2)) then
1293  write(ncpout2_cll,*) ' Further output of Critical Points for TBten_cll suppressed'
1294  write(ncpout2_cll,*)
1295  endif
1296 #endif
1297  end if
1298  end if
1299  end if
1300 

◆ bten_args_list_cll()

subroutine collier_tensors::bten_cll::bten_args_list_cll ( double complex, dimension(:), intent(out)  TB,
double complex, dimension(:), intent(out)  TBuv,
double complex, dimension(0:3), intent(in)  p1vec,
double complex, intent(in)  p10,
double complex, intent(in)  m02,
double complex, intent(in)  m12,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  TBerr 
)

Definition at line 1313 of file collier_tensors.F90.

1313 
1314  integer, intent(in) :: rmax
1315  double complex, intent(in) :: p1vec(0:3)
1316  double complex, intent(in) :: p10,m02,m12
1317  double complex, intent(out) :: TB(:), TBuv(:)
1318  double precision, intent(out), optional :: TBerr(0:rmax)
1319  integer :: r,i
1320  logical :: eflag
1321 
1322  if (2.gt.nmax_cll) then
1323  call seterrflag_cll(-10)
1324  call errout_cll('Bten_cll','Nmax_cll smaller 2',eflag,.true.)
1325  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
1326  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 2'
1327  call propagateerrflag_cll
1328  return
1329  end if
1330  if (rmax.gt.rmax_cll) then
1331  call seterrflag_cll(-10)
1332  call errout_cll('Bten_cll','argument rmax larger than rmax_cll',eflag,.true.)
1333  write(nerrout_cll,*) 'rmax =',rmax
1334  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
1335  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
1336  call propagateerrflag_cll
1337  return
1338  end if
1339 
1340  call bten_args_list_checked_cll(tb,tbuv,p1vec,p10,m02,m12,rmax,tberr)
1341 

◆ bten_list_cll()

subroutine collier_tensors::bten_cll::bten_list_cll ( double complex, dimension(:), intent(out)  TB,
double complex, dimension(:), intent(out)  TBuv,
double complex, dimension(0:3,1), intent(in)  MomVec,
double complex, dimension(1), intent(in)  MomInv,
double complex, dimension(0:1), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  TBerr 
)

Definition at line 962 of file collier_tensors.F90.

962 
963  integer, intent(in) :: rmax
964  double complex, intent(in) :: MomVec(0:3,1), MomInv(1), masses2(0:1)
965  double complex, intent(out) :: TB(:), TBuv(:)
966  double precision, intent(out), optional :: TBerr(0:rmax)
967  integer :: r,i
968  logical :: eflag
969 
970  if (2.gt.nmax_cll) then
971  call seterrflag_cll(-10)
972  call errout_cll('Bten_cll','Nmax_cll smaller 2',eflag,.true.)
973  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
974  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 2'
975  call propagateerrflag_cll
976  return
977  end if
978  if (rmax.gt.rmax_cll) then
979  call seterrflag_cll(-10)
980  call errout_cll('Bten_cll','argument rmax larger than rmax_cll',eflag,.true.)
981  write(nerrout_cll,*) 'rmax =',rmax
982  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
983  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
984  call propagateerrflag_cll
985  return
986  end if
987 
988  call bten_list_checked_cll(tb,tbuv,momvec,mominv,masses2,rmax,tberr)
989 

◆ bten_main_cll()

subroutine collier_tensors::bten_cll::bten_main_cll ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TB,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TBuv,
double complex, dimension(0:3,1), intent(in)  MomVec,
double complex, dimension(1), intent(in)  MomInv,
double complex, dimension(0:1), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  TBerr 
)

Definition at line 791 of file collier_tensors.F90.

791 
792  integer, intent(in) :: rmax
793  double complex, intent(in) :: MomVec(0:3,1), MomInv(1), masses2(0:1)
794  double complex, intent(out) :: TB(0:rmax,0:rmax,0:rmax,0:rmax)
795  double complex, intent(out) :: TBuv(0:rmax,0:rmax,0:rmax,0:rmax)
796  double precision, intent(out), optional :: TBerr(0:rmax)
797  double complex :: TB2(0:rmax,0:rmax,0:rmax,0:rmax), TBuv2(0:rmax,0:rmax,0:rmax,0:rmax)
798  double complex :: CB(0:rmax/2,0:rmax), CBuv(0:rmax/2,0:rmax)
799  double precision :: CBerr(0:rmax), TBerr_aux(0:rmax), TBerr_aux2(0:rmax)
800  double complex :: args(7)
801  double precision :: TBdiff(0:rmax),norm(0:rmax),norm_coli,norm_dd,TBacc(0:rmax)
802  integer :: r,n0,n1,n2,n3
803  logical :: eflag
804 
805  if (2.gt.nmax_cll) then
806  call seterrflag_cll(-10)
807  call errout_cll('Bten_cll','Nmax_cll smaller 2',eflag,.true.)
808  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
809  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 2'
810  call propagateerrflag_cll
811  return
812  end if
813  if (rmax.gt.rmax_cll) then
814  call seterrflag_cll(-10)
815  call errout_cll('Bten_cll','argument rmax larger than rmax_cll',eflag,.true.)
816  write(nerrout_cll,*) 'rmax =',rmax
817  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
818  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
819  call propagateerrflag_cll
820  return
821  end if
822 
823  ! set ID of master call
824  args(1:4) = momvec(0:,1)
825  args(5) = mominv(1)
826  args(6:7) = masses2(0:)
827  call setmasterfname_cll('Bten_cll')
828  call setmastern_cll(2)
829  call setmasterr_cll(rmax)
830  call setmasterargs_cll(7,args)
831 
832  call settencache_cll(tenred_cll-1)
833 
834  if (mode_cll.eq.3) then
835  ! calculate tensor with coefficients from COLI
836  mode_cll = 1
837  call b_main_cll(cb,cbuv,mominv(1),masses2(0),masses2(1),rmax,cberr,0)
838  call calctensorb(tb,tbuv,tberr_aux,cb,cbuv,cberr,momvec(0:,1),rmax)
839 
840  ! calculate tensor with coefficients from DD
841  mode_cll = 2
842  call b_main_cll(cb,cbuv,mominv(1),masses2(0),masses2(1),rmax,cberr,0)
843  call calctensorb(tb2,tbuv2,tberr_aux2,cb,cbuv,cberr,momvec(0:,1),rmax)
844 
845  ! comparison --> take better result
846  mode_cll = 3
847  do r=0,rmax
848  norm_coli=0d0
849  norm_dd=0d0
850  do n0=0,r
851  do n1=0,r-n0
852  do n2=0,r-n0-n1
853  n3=r-n0-n1-n2
854  norm_coli = max(norm_coli,abs(tb(n0,n1,n2,n3)))
855  norm_dd = max(norm_dd,abs(tb2(n0,n1,n2,n3)))
856  end do
857  end do
858  end do
859  if (norm_coli.eq.0d0) then
860  norm_coli = max(abs(mominv(1)),maxval(abs(masses2(0:1))))
861  if(norm_coli.ne.0d0) then
862  norm_coli=norm_coli**(real(r)/2)
863  else
864  norm_coli=muir2_cll**(real(r)/2)
865  end if
866  end if
867  if (norm_dd.eq.0d0) then
868  norm_dd = max(abs(mominv(1)),maxval(abs(masses2(0:1))))
869  if(norm_dd.ne.0d0) then
870  norm_dd=norm_dd**(real(r)/2)
871  else
872  norm_dd=muir2_cll**(real(r)/2)
873  end if
874  end if
875  norm(r) = min(norm_coli,norm_dd)
876  end do
877 
878  call checktensors_cll(tb,tb2,momvec,mominv,masses2,norm,2,rmax,tbdiff)
879 
880  if (tberr_aux(rmax).lt.tberr_aux2(rmax)) then
881  if (present(tberr)) tberr = max(tberr_aux,tbdiff*norm)
882  do r=0,rmax
883  tbacc(r) = max(tberr_aux(r)/norm(r),tbdiff(r))
884  end do
885  if (monitoring) pointscntbten_coli = pointscntbten_coli + 1
886  else
887  tb = tb2
888  tbuv = tbuv2
889  if (present(tberr)) tberr = max(tberr_aux2,tbdiff*norm)
890  do r=0,rmax
891  tbacc(r) = max(tberr_aux2(r)/norm(r),tbdiff(r))
892  end do
893  if (monitoring) pointscntbten_dd = pointscntbten_dd + 1
894  end if
895 
896  else
897  call b_main_cll(cb,cbuv,mominv(1),masses2(0),masses2(1),rmax,cberr,0)
898  call calctensorb(tb,tbuv,tberr_aux,cb,cbuv,cberr,momvec(0:,1),rmax)
899  if (present(tberr)) tberr = tberr_aux
900  norm = 0d0
901  do r=0,rmax
902  do n0=0,r
903  do n1=0,r-n0
904  do n2=0,r-n0-n1
905  n3=r-n0-n1-n2
906  norm(r) = max(norm(r),abs(tb(n0,n1,n2,n3)))
907  end do
908  end do
909  end do
910  if (norm(r).eq.0d0) then
911  norm(r) = max(abs(mominv(1)),maxval(abs(masses2(0:1))))
912  if(norm(r).ne.0d0) then
913  norm(r)=norm(r)**(real(r)/2)
914  else
915  norm(r)=muir2_cll**(real(r)/2)
916  end if
917  end if
918  tbacc(r) = tberr_aux(r)/norm(r)
919  end do
920 
921  end if
922 
923  call propagateaccflag_cll(tbacc,rmax)
924  call propagateerrflag_cll
925 
926  if (monitoring) then
927  pointscntbten_cll = pointscntbten_cll + 1
928 
929  if(maxval(tbacc).gt.reqacc_cll) accpointscntbten_cll = accpointscntbten_cll + 1
930 
931  if(maxval(tbacc).gt.critacc_cll) then
932  critpointscntbten_cll = critpointscntbten_cll + 1
933  if ( critpointscntbten_cll.le.noutcritpointsmax_cll(2) ) then
934  call critpointsout_cll('TBten_cll',0,maxval(tbacc),critpointscntbten_cll)
935  if( critpointscntbten_cll.eq.noutcritpointsmax_cll(2)) then
936  write(ncpout_cll,*) ' Further output of Critical Points for TBten_cll suppressed'
937  write(ncpout_cll,*)
938  endif
939 #ifdef CritPoints2
940  call critpointsout2_cll('TBten_cll',0,maxval(tbacc),critpointscntbten_cll)
941  if( critpointscntbten_cll.eq.noutcritpointsmax_cll(2)) then
942  write(ncpout2_cll,*) ' Further output of Critical Points for TBten_cll suppressed'
943  write(ncpout2_cll,*)
944  endif
945 #endif
946  end if
947  end if
948  end if
949 

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