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.
Public Member Functions | List of all members
collier_coefs::d_cll Interface Reference

Public Member Functions

subroutine d_main_cll (D, Duv, p10, p21, p32, p30, p20, p31, m02, m12, m22, m32, rmax, Derr, id_in, Derr2)
 
subroutine d_arrays_cll (D, Duv, MomInv, masses2, rmax, Derr, Derr2)
 
subroutine d_list_cll (D, Duv, p10, p21, p32, p30, p20, p31, m02, m12, m22, m32, rmax, Derr, Derr2)
 
subroutine d_arrays_list_cll (D, Duv, MomInv, masses2, rmax, Derr, Derr2)
 

Detailed Description

Definition at line 52 of file collier_coefs.F90.

Member Function/Subroutine Documentation

◆ d_arrays_cll()

subroutine collier_coefs::d_cll::d_arrays_cll ( double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax), intent(out)  D,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax), intent(out)  Duv,
double complex, dimension(6), intent(in)  MomInv,
double complex, dimension(0:3), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Derr,
double precision, dimension(0:rmax), intent(out), optional  Derr2 
)

Definition at line 1545 of file collier_coefs.F90.

1545 
1546  integer, intent(in) :: rmax
1547  double complex, intent(in) :: MomInv(6), masses2(0:3)
1548  double complex, intent(out) :: D(0:rmax/2,0:rmax,0:rmax,0:rmax)
1549  double complex, intent(out) :: Duv(0:rmax/2,0:rmax,0:rmax,0:rmax)
1550  double precision, optional, intent(out) :: Derr(0:rmax),Derr2(0:rmax)
1551  double precision :: Derraux(0:rmax),Derr2aux(0:rmax)
1552 
1553  logical :: eflag
1554 
1555  if (4.gt.nmax_cll) then
1556  call seterrflag_cll(-10)
1557  call errout_cll('D_cll','Nmax_cll smaller 4',eflag,.true.)
1558  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
1559  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 4'
1560  call propagateerrflag_cll
1561  return
1562  end if
1563  if (rmax.gt.rmax_cll) then
1564  call seterrflag_cll(-10)
1565  call errout_cll('D_cll','argument rmax larger than rmax_cll',eflag,.true.)
1566  write(nerrout_cll,*) 'rmax =',rmax
1567  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
1568  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
1569  call propagateerrflag_cll
1570  return
1571  end if
1572 
1573  if (present(derr)) then
1574  if (present(derr2)) then
1575  call d_main_cll(d,duv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
1576  masses2(0),masses2(1),masses2(2),masses2(3),rmax,derr,derr2=derr2)
1577  else
1578  call d_main_cll(d,duv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
1579  masses2(0),masses2(1),masses2(2),masses2(3),rmax,derr)
1580  end if
1581  else
1582  if (present(derr2)) then
1583  call d_main_cll(d,duv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
1584  masses2(0),masses2(1),masses2(2),masses2(3),rmax,derraux,derr2=derr2)
1585  else
1586  call d_main_cll(d,duv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
1587  masses2(0),masses2(1),masses2(2),masses2(3),rmax,derraux)
1588  end if
1589  end if
1590 

◆ d_arrays_list_cll()

subroutine collier_coefs::d_cll::d_arrays_list_cll ( double complex, dimension(:), intent(out)  D,
double complex, dimension(:), intent(out)  Duv,
double complex, dimension(6), intent(in)  MomInv,
double complex, dimension(0:3), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Derr,
double precision, dimension(0:rmax), intent(out), optional  Derr2 
)

Definition at line 1691 of file collier_coefs.F90.

1691 
1692  integer, intent(in) :: rmax
1693  double complex, intent(in) :: MomInv(6), masses2(0:3)
1694  double complex, intent(out) :: D(:),Duv(:)
1695  double precision, optional, intent(out) :: Derr(0:rmax),Derr2(0:rmax)
1696  logical :: eflag
1697 
1698  if (4.gt.nmax_cll) then
1699  call seterrflag_cll(-10)
1700  call errout_cll('D_cll','Nmax_cll smaller 4',eflag,.true.)
1701  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
1702  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 4'
1703  call propagateerrflag_cll
1704  return
1705  end if
1706  if (rmax.gt.rmax_cll) then
1707  call seterrflag_cll(-10)
1708  call errout_cll('D_cll','argument rmax larger than rmax_cll',eflag,.true.)
1709  write(nerrout_cll,*) 'rmax =',rmax
1710  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
1711  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
1712  call propagateerrflag_cll
1713  return
1714  end if
1715 
1716  call d_arrays_list_checked_cll(d,duv,mominv,masses2,rmax,derr,derr2)
1717 

◆ d_list_cll()

subroutine collier_coefs::d_cll::d_list_cll ( double complex, dimension(:), intent(out)  D,
double complex, dimension(:), intent(out)  Duv,
double complex, intent(in)  p10,
double complex, intent(in)  p21,
double complex, intent(in)  p32,
double complex, intent(in)  p30,
double complex, intent(in)  p20,
double complex, intent(in)  p31,
double complex, intent(in)  m02,
double complex, intent(in)  m12,
double complex, intent(in)  m22,
double complex, intent(in)  m32,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Derr,
double precision, dimension(0:rmax), intent(out), optional  Derr2 
)

Definition at line 1603 of file collier_coefs.F90.

1603 
1604  integer, intent(in) :: rmax
1605  double complex, intent(in) :: p10,p21,p32,p30,p20,p31,m02,m12,m22,m32
1606  double complex, intent(out) :: D(:),Duv(:)
1607  double precision, optional, intent(out) :: Derr(0:rmax),Derr2(0:rmax)
1608  logical :: eflag
1609 
1610  if (4.gt.nmax_cll) then
1611  call seterrflag_cll(-10)
1612  call errout_cll('D_cll','Nmax_cll smaller 4',eflag,.true.)
1613  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
1614  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 4'
1615  call propagateerrflag_cll
1616  return
1617  end if
1618  if (rmax.gt.rmax_cll) then
1619  call seterrflag_cll(-10)
1620  call errout_cll('D_cll','argument rmax larger than rmax_cll',eflag,.true.)
1621  write(nerrout_cll,*) 'rmax =',rmax
1622  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
1623  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
1624  call propagateerrflag_cll
1625  return
1626  end if
1627 
1628  call d_list_checked_cll(d,duv,p10,p21,p32,p30,p20,p31, &
1629  m02,m12,m22,m32,rmax,derr,derr2)
1630 

◆ d_main_cll()

subroutine collier_coefs::d_cll::d_main_cll ( double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax), intent(out)  D,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax), intent(out)  Duv,
double complex, intent(in)  p10,
double complex, intent(in)  p21,
double complex, intent(in)  p32,
double complex, intent(in)  p30,
double complex, intent(in)  p20,
double complex, intent(in)  p31,
double complex, intent(in)  m02,
double complex, intent(in)  m12,
double complex, intent(in)  m22,
double complex, intent(in)  m32,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Derr,
integer, intent(in), optional  id_in,
double precision, dimension(0:rmax), intent(out), optional  Derr2 
)

Definition at line 1213 of file collier_coefs.F90.

1213 
1214  integer, intent(in) :: rmax
1215  double complex, intent(in) :: p10,p21,p32,p30,p20,p31,m02,m12,m22,m32
1216  double precision :: q10,q21,q32,q30,q20,q31
1217  double complex :: mm02,mm12,mm22,mm32
1218  double complex, intent(out) :: D(0:rmax/2,0:rmax,0:rmax,0:rmax)
1219  double complex, intent(out) :: Duv(0:rmax/2,0:rmax,0:rmax,0:rmax)
1220  double precision, optional, intent(out) :: Derr(0:rmax),Derr2(0:rmax)
1221  integer, optional, intent(in) :: id_in
1222  double complex :: D2uv(0:rmax/2,0:rmax,0:rmax,0:rmax)
1223  double complex :: D2(0:rmax/2,0:rmax,0:rmax,0:rmax)
1224  double complex :: Dcoliuv(0:rmax,0:rmax,0:rmax,0:rmax)
1225  double complex :: Dcoli(0:rmax,0:rmax,0:rmax,0:rmax)
1226  double complex :: Ddduv(0:rmax,0:rmax,0:rmax,0:rmax)
1227  double complex :: Ddd(0:rmax,0:rmax,0:rmax,0:rmax)
1228  double precision :: Derraux(0:rmax),Derr2aux(0:rmax),Ddiff(0:rmax)
1229  double complex :: elimcminf2
1230  double complex :: args(10)
1231  integer :: n0,rank,errflag,id
1232  double precision :: accrelDD(0:rmax_DD),accabsDD(0:rmax_DD),Dacc(0:rmax),norm,norm_coli,norm_dd,Dacc2(0:rmax)
1233  double precision :: accrel2DD(0:rmax_DD),accabs2DD(0:rmax_DD)
1234  integer :: accflagDD,errflagDD,NDD,rankDD
1235  logical :: mflag,eflag
1236  integer :: r,n1,n2,n3
1237 
1238  if (4.gt.nmax_cll) then
1239  call seterrflag_cll(-10)
1240  call errout_cll('D_cll','Nmax_cll smaller 4',eflag,.true.)
1241  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
1242  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 4'
1243  call propagateerrflag_cll
1244  return
1245  end if
1246  if (rmax.gt.rmax_cll) then
1247  call seterrflag_cll(-10)
1248  call errout_cll('D_cll','argument rmax larger than rmax_cll',eflag,.true.)
1249  write(nerrout_cll,*) 'rmax =',rmax
1250  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
1251  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
1252  call propagateerrflag_cll
1253  return
1254  end if
1255 
1256  mflag=.true.
1257  if (present(id_in)) then
1258  mflag=.false.
1259  id = id_in
1260  else
1261  id = 0
1262  end if
1263 
1264  if (mflag) then
1265  ! set ID of master call
1266  args(1) = p10
1267  args(2) = p21
1268  args(3) = p32
1269  args(4) = p30
1270  args(5) = p20
1271  args(6) = p31
1272  args(7) = m02
1273  args(8) = m12
1274  args(9) = m22
1275  args(10) = m32
1276  call setmasterfname_cll('D_cll')
1277  call setmastern_cll(4)
1278  call setmasterr_cll(rmax)
1279  call setmasterargs_cll(10,args)
1280 
1281  call settencache_cll(never_tenred_cll)
1282  end if
1283 
1284 
1285  select case (mode_cll)
1286 
1287  case (1)
1288  ! calculate loop integral using
1289  ! COLI implementation by AD/LH
1290 
1291  call calcd(dcoli,dcoliuv,p10,p21,p32,p30,p20,p31, &
1292  m02,m12,m22,m32,rmax,id,derraux,derr2aux)
1293 
1294  norm = abs(dcoli(0,0,0,0))
1295  do r=1,rmax
1296  do n1=0,r
1297  do n2=0,r-n1
1298  n3=r-n1-n2
1299  norm = max(norm,abs(dcoli(0,n1,n2,n3)))
1300  end do
1301  end do
1302  end do
1303  if (norm.eq.0d0) then
1304  norm = max(abs(p10),abs(p21),abs(p32),abs(p30),abs(p20),abs(p31), &
1305  abs(m02),abs(m12),abs(m22),abs(m32))
1306  if(norm.ne.0d0) then
1307  norm=1d0/norm**2
1308  else
1309  norm=1d0/muir2_cll**2
1310  end if
1311  end if
1312  if (norm.ne.0d0) then
1313  dacc = derraux/norm
1314  dacc2 = derr2aux/norm
1315  else
1316  dacc = 0d0
1317  dacc2 = 0d0
1318  end if
1319 
1320  if (present(derr)) derr = derraux
1321  if (present(derr2)) derr2 = derr2aux
1322 
1323  if (mflag) call propagateaccflag_cll(dacc,rmax)
1324 
1325  d(0:rmax/2,0:rmax,0:rmax,0:rmax) = dcoli(0:rmax/2,0:rmax,0:rmax,0:rmax)
1326  duv(0:rmax/2,0:rmax,0:rmax,0:rmax) = dcoliuv(0:rmax/2,0:rmax,0:rmax,0:rmax)
1327 
1328 
1329  case (2)
1330  ! calculate loop integral using
1331  ! DD implementation by SD
1332 
1333  id=0
1334 
1335  ! replace small masses by DD-identifiers
1336  q10 = dreal(getminf2dd_cll(p10))
1337  q21 = dreal(getminf2dd_cll(p21))
1338  q32 = dreal(getminf2dd_cll(p32))
1339  q30 = dreal(getminf2dd_cll(p30))
1340  q20 = dreal(getminf2dd_cll(p20))
1341  q31 = dreal(getminf2dd_cll(p31))
1342  mm02 = getminf2dd_cll(m02)
1343  mm12 = getminf2dd_cll(m12)
1344  mm22 = getminf2dd_cll(m22)
1345  mm32 = getminf2dd_cll(m32)
1346 
1347  rank = rmax
1348 ! write(*,*) rank
1349 ! write(*,*) q10,q21
1350 ! write(*,*) q32,q30
1351 ! write(*,*) q20,q31
1352 ! write(*,*) mm02,mm12
1353 ! write(*,*) mm22,mm32
1354  call d_dd(ddd,ddduv,q10,q21,q32,q30,q20,q31, &
1355  mm02,mm12,mm22,mm32,rank,id)
1356  d(0:rank/2,0:rank,0:rank,0:rank) = ddd(0:rank/2,0:rank,0:rank,0:rank)
1357  duv(0:rank/2,0:rank,0:rank,0:rank) = ddduv(0:rank/2,0:rank,0:rank,0:rank)
1358 
1359  call ddgetacc(accreldd,accabsdd,accrel2dd,accabs2dd,ndd,rankdd,accflagdd,errflagdd,id)
1360  if (present(derr)) derr(0:rmax) = accabsdd(0:rmax)
1361  if (present(derr2)) derr2(0:rmax) = accabs2dd(0:rmax)
1362 
1363  norm = abs(d(0,0,0,0))
1364  do r=1,rmax
1365  do n1=0,r
1366  do n2=0,r-n1
1367  n3=r-n1-n2
1368  norm = max(norm,abs(d(0,n1,n2,n3)))
1369  end do
1370  end do
1371  end do
1372  if (norm.eq.0d0) then
1373  norm = max(abs(p10),abs(p21),abs(p32),abs(p30),abs(p20),abs(p31), &
1374  abs(m02),abs(m12),abs(m22),abs(m32))
1375  if(norm.ne.0d0) then
1376  norm=1d0/norm**2
1377  else
1378  norm=1d0/muir2_cll**2
1379  end if
1380  end if
1381  if (norm.ne.0d0) then
1382  dacc = accabsdd(0:rmax)/norm
1383  dacc2 = accabs2dd(0:rmax)/norm
1384  else
1385  dacc = 0d0
1386  dacc2 = 0d0
1387  end if
1388  if (mflag) call propagateaccflag_cll(dacc,rmax)
1389 
1390 
1391  case (3)
1392  ! cross-check mode
1393  ! compare results for loop integral
1394  ! from COLI implementation by AD/LH and
1395  ! from DD implementation by SD
1396 
1397  ! calculate loop integral using COLI
1398  call calcd(dcoli,dcoliuv,p10,p21,p32,p30,p20,p31, &
1399  m02,m12,m22,m32,rmax,id,derraux,derr2aux)
1400 
1401  d(0:rmax/2,0:rmax,0:rmax,0:rmax) = dcoli(0:rmax/2,0:rmax,0:rmax,0:rmax)
1402  duv(0:rmax/2,0:rmax,0:rmax,0:rmax) = dcoliuv(0:rmax/2,0:rmax,0:rmax,0:rmax)
1403 
1404 
1405  ! calculate loop integral using DD
1406 
1407  id=0
1408 
1409  ! replace small masses by DD-identifiers
1410  q10 = dreal(getminf2dd_cll(p10))
1411  q21 = dreal(getminf2dd_cll(p21))
1412  q32 = dreal(getminf2dd_cll(p32))
1413  q30 = dreal(getminf2dd_cll(p30))
1414  q20 = dreal(getminf2dd_cll(p20))
1415  q31 = dreal(getminf2dd_cll(p31))
1416  mm02 = getminf2dd_cll(m02)
1417  mm12 = getminf2dd_cll(m12)
1418  mm22 = getminf2dd_cll(m22)
1419  mm32 = getminf2dd_cll(m32)
1420 
1421  rank = rmax
1422  call d_dd(ddd,ddduv,q10,q21,q32,q30,q20,q31, &
1423  mm02,mm12,mm22,mm32,rank,id)
1424  do n0=0,rank/2
1425  d2(n0,0:rank,0:rank,0:rank) = ddd(n0,0:rank,0:rank,0:rank)
1426  d2uv(n0,0:rank,0:rank,0:rank) = ddduv(n0,0:rank,0:rank,0:rank)
1427  end do
1428  call ddgetacc(accreldd,accabsdd,accrel2dd,accabs2dd,ndd,rankdd,accflagdd,errflagdd,id)
1429 
1430  norm_coli = abs(d(0,0,0,0))
1431  norm_dd = abs(d2(0,0,0,0))
1432  do r=1,rmax
1433  do n1=0,r
1434  do n2=0,r-n1
1435  n3=r-n1-n2
1436  norm_coli = max(norm_coli,abs(d(0,n1,n2,n3)))
1437  norm_dd = max(norm_dd,abs(d2(0,n1,n2,n3)))
1438  end do
1439  end do
1440  end do
1441  if (norm_coli.eq.0d0) then
1442  norm_coli = max(abs(p10),abs(p21),abs(p32),abs(p30),abs(p20),abs(p31), &
1443  abs(m02),abs(m12),abs(m22),abs(m32))
1444  if(norm_coli.ne.0d0) then
1445  norm_coli=1d0/norm_coli**2
1446  else
1447  norm_coli=1d0/muir2_cll**2
1448  end if
1449  end if
1450  if (norm_dd.eq.0d0) then
1451  norm_dd = max(abs(p10),abs(p21),abs(p32),abs(p30),abs(p20),abs(p31), &
1452  abs(m02),abs(m12),abs(m22),abs(m32))
1453  if(norm_dd.ne.0d0) then
1454  norm_dd=1d0/norm_dd**2
1455  else
1456  norm_dd=1d0/muir2_cll**2
1457  end if
1458  end if
1459  norm = min(norm_coli,norm_dd)
1460 
1461  ! cross-check
1462  call checkcoefsd_cll(d,d2,p10,p21,p32,p30,p20,p31, &
1463  m02,m12,m22,m32,rmax,norm,ddiff)
1464 
1465 
1466  if (derraux(rmax).lt.accabsdd(rmax)) then
1467  if (present(derr)) derr = max(derraux,ddiff)
1468  if (present(derr2)) derr2 = derr2aux
1469  if (norm.ne.0d0) then
1470  dacc = max(derraux/norm_coli,ddiff/norm)
1471  dacc2 = derr2aux/norm_coli
1472  else
1473  dacc = ddiff
1474  dacc2 = 0d0
1475  end if
1476  if (monitoring) pointscntd_coli = pointscntd_coli + 1
1477  else
1478  d = d2
1479  duv = d2uv
1480  if (present(derr)) derr = max(accabsdd(0:rmax),ddiff)
1481  if (present(derr2)) derr2 = accabs2dd(0:rmax)
1482  if (norm.ne.0d0) then
1483  dacc = max(accabsdd(0:rmax)/norm_dd,ddiff/norm)
1484  dacc2 = accabs2dd(0:rmax)/norm_dd
1485  else
1486  dacc = ddiff
1487  dacc2 = 0d0
1488  end if
1489  if (monitoring) pointscntd_dd = pointscntd_dd + 1
1490  end if
1491 
1492  if (mflag) call propagateaccflag_cll(dacc,rmax)
1493 
1494  end select
1495 
1496  if (mflag) call propagateerrflag_cll
1497 
1498  if (monitoring) then
1499  pointscntd_cll = pointscntd_cll + 1
1500 
1501  if(maxval(dacc).gt.reqacc_cll) accpointscntd_cll = accpointscntd_cll + 1
1502 
1503  if(maxval(dacc).gt.critacc_cll) then
1504  critpointscntd_cll = critpointscntd_cll + 1
1505  if ( critpointscntd_cll.le.noutcritpointsmax_cll(4) ) then
1506  call critpointsout_cll('D_cll',0,maxval(dacc), critpointscntd_cll)
1507  if( critpointscntd_cll.eq.noutcritpointsmax_cll(4)) then
1508  write(ncpout_cll,*) ' Further output of Critical Points for D_cll suppressed '
1509  write(ncpout_cll,*)
1510  endif
1511  end if
1512  end if
1513 
1514 
1515 #ifdef CritPoints2
1516  if(maxval(dacc2).gt.reqacc_cll) accpointscntd2_cll = accpointscntd2_cll + 1
1517 
1518  if(maxval(dacc2).gt.critacc_cll) then
1519  critpointscntd2_cll = critpointscntd2_cll + 1
1520  if ( critpointscntd2_cll.le.noutcritpointsmax_cll(4) ) then
1521  call critpointsout2_cll('D_cll',0,maxval(dacc2), critpointscntd2_cll)
1522  if( critpointscntd2_cll.eq.noutcritpointsmax_cll(4)) then
1523  write(ncpout2_cll,*) ' Further output of Critical Points for D_cll suppressed '
1524  write(ncpout2_cll,*)
1525  endif
1526  end if
1527  end if
1528 #endif
1529 
1530  end if
1531 
1532 

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