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

Public Member Functions

subroutine cten_main_cll (TC, TCuv, MomVec, MomInv, masses2, rmax, TCerr)
 
subroutine cten_list_cll (TC, TCuv, MomVec, MomInv, masses2, rmax, TCerr)
 
subroutine cten_args_cll (TC, TCuv, p1vec, p2vec, p10, p21, p20, m02, m12, m22, rmax, TCerr)
 
subroutine cten_args_list_cll (TC, TCuv, p1vec, p2vec, p10, p21, p20, m02, m12, m22, rmax, TCerr)
 

Detailed Description

Definition at line 56 of file collier_tensors.F90.

Member Function/Subroutine Documentation

◆ cten_args_cll()

subroutine collier_tensors::cten_cll::cten_args_cll ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TC,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TCuv,
double complex, dimension(0:3), intent(in)  p1vec,
double complex, dimension(0:3), intent(in)  p2vec,
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  TCerr 
)

Definition at line 1846 of file collier_tensors.F90.

1846 
1847  integer, intent(in) :: rmax
1848  double complex, intent(in) :: p1vec(0:3), p2vec(0:3)
1849  double complex, intent(in) :: p10,p21,p20,m02,m12,m22
1850  double complex, intent(out) :: TC(0:rmax,0:rmax,0:rmax,0:rmax)
1851  double complex, intent(out) :: TCuv(0:rmax,0:rmax,0:rmax,0:rmax)
1852  double precision, intent(out), optional :: TCerr(0:rmax)
1853  double complex :: TC2(0:rmax,0:rmax,0:rmax,0:rmax), TCuv2(0:rmax,0:rmax,0:rmax,0:rmax)
1854  double complex :: MomVec(0:3,2), MomInv(3), masses2(0:2)
1855  double complex :: CC(0:rmax/2,0:rmax,0:rmax), CCuv(0:rmax/2,0:rmax,0:rmax)
1856  double precision :: CCerr(0:rmax), TCerr_aux(0:rmax), TCerr_aux2(0:rmax)
1857  double complex :: args(14)
1858  double precision :: TCdiff(0:rmax),norm(0:rmax),norm_coli,norm_dd,TCacc(0:rmax)
1859  integer :: r,n0,n1,n2,n3
1860  logical :: eflag
1861 
1862  if (3.gt.nmax_cll) then
1863  call seterrflag_cll(-10)
1864  call errout_cll('Cten_cll','Nmax_cll smaller 3',eflag,.true.)
1865  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
1866  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 3'
1867  call propagateerrflag_cll
1868  return
1869  end if
1870  if (rmax.gt.rmax_cll) then
1871  call seterrflag_cll(-10)
1872  call errout_cll('Cten_cll','argument rmax larger than rmax_cll',eflag,.true.)
1873  write(nerrout_cll,*) 'rmax =',rmax
1874  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
1875  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
1876  call propagateerrflag_cll
1877  return
1878  end if
1879 
1880  momvec(0:,1) = p1vec
1881  momvec(0:,2) = p2vec
1882  mominv(1) = p10
1883  mominv(2) = p21
1884  mominv(3) = p20
1885  masses2(0) = m02
1886  masses2(1) = m12
1887  masses2(2) = m22
1888 
1889  ! set ID of master call
1890  args(1:4) = momvec(0:,1)
1891  args(5:8) = momvec(0:,2)
1892  args(9:11) = mominv
1893  args(12:14) = masses2(0:)
1894  call setmasterfname_cll('Cten_cll')
1895  call setmastern_cll(3)
1896  call setmasterr_cll(rmax)
1897  call setmasterargs_cll(14,args)
1898 
1899  call settencache_cll(tenred_cll-1)
1900 
1901 
1902  if (mode_cll.eq.3) then
1903  ! calculate tensor with coefficients from COLI
1904  mode_cll = 1
1905  call c_main_cll(cc,ccuv,mominv(1),mominv(2),mominv(3), &
1906  masses2(0),masses2(1),masses2(2),rmax,cerr2=ccerr,id_in=0)
1907  call calctensorc(tc,tcuv,tcerr_aux,cc,ccuv,ccerr,momvec,rmax)
1908 
1909  ! calculate tensor with coefficients from DD
1910  mode_cll = 2
1911  call c_main_cll(cc,ccuv,mominv(1),mominv(2),mominv(3), &
1912  masses2(0),masses2(1),masses2(2),rmax,cerr2=ccerr,id_in=0)
1913  call calctensorc(tc2,tcuv2,tcerr_aux2,cc,ccuv,ccerr,momvec,rmax)
1914 
1915  ! comparison --> take better result
1916  mode_cll = 3
1917  do r=0,rmax
1918  norm_coli=0d0
1919  norm_dd=0d0
1920  do n0=0,r
1921  do n1=0,r-n0
1922  do n2=0,r-n0-n1
1923  n3=r-n0-n1-n2
1924  norm_coli = max(norm_coli,abs(tc(n0,n1,n2,n3)))
1925  norm_dd = max(norm_dd,abs(tc2(n0,n1,n2,n3)))
1926  end do
1927  end do
1928  end do
1929  if (norm_coli.eq.0d0) then
1930  norm_coli = max(maxval(abs(mominv(1:3))),maxval(abs(masses2(0:2))))
1931  if(norm_coli.ne.0d0) then
1932  norm_coli=1d0/norm_coli**(1-real(r)/2)
1933  else
1934  norm_coli=1d0/muir2_cll**(1-real(r)/2)
1935  end if
1936  end if
1937  if (norm_dd.eq.0d0) then
1938  norm_dd = max(maxval(abs(mominv(1:3))),maxval(abs(masses2(0:2))))
1939  if(norm_dd.ne.0d0) then
1940  norm_dd=1d0/norm_dd**(1-real(r)/2)
1941  else
1942  norm_dd=1d0/muir2_cll**(1-real(r)/2)
1943  end if
1944  end if
1945  norm(r) = min(norm_coli,norm_dd)
1946  end do
1947 
1948  call checktensors_cll(tc,tc2,momvec,mominv,masses2,norm,3,rmax,tcdiff)
1949 
1950  if (tcerr_aux(rmax).lt.tcerr_aux2(rmax)) then
1951  if (present(tcerr)) tcerr = max(tcerr_aux,tcdiff*norm)
1952  do r=0,rmax
1953  tcacc(r) = max(tcerr_aux(r)/norm(r),tcdiff(r))
1954  end do
1955  if (monitoring) pointscntcten_coli = pointscntcten_coli + 1
1956  else
1957  tc = tc2
1958  tcuv = tcuv2
1959  if (present(tcerr)) tcerr = max(tcerr_aux2,tcdiff*norm)
1960  do r=0,rmax
1961  tcacc(r) = max(tcerr_aux2(r)/norm(r),tcdiff(r))
1962  end do
1963  if (monitoring) pointscntcten_dd = pointscntcten_dd + 1
1964  end if
1965 
1966  else
1967  call c_main_cll(cc,ccuv,mominv(1),mominv(2),mominv(3), &
1968  masses2(0),masses2(1),masses2(2),rmax,cerr2=ccerr,id_in=0)
1969  call calctensorc(tc,tcuv,tcerr_aux,cc,ccuv,ccerr,momvec,rmax)
1970  if (present(tcerr)) tcerr = tcerr_aux
1971  norm=0d0
1972  do r=0,rmax
1973  do n0=0,r
1974  do n1=0,r-n0
1975  do n2=0,r-n0-n1
1976  n3=r-n0-n1-n2
1977  norm(r) = max(norm(r),abs(tc(n0,n1,n2,n3)))
1978  end do
1979  end do
1980  end do
1981  if (norm(r).eq.0d0) then
1982  norm(r) = max(maxval(abs(mominv(1:3))),maxval(abs(masses2(0:2))))
1983  if(norm(r).ne.0d0) then
1984  norm(r)=1d0/norm(r)**(1-real(r)/2)
1985  else
1986  norm(r)=1d0/muir2_cll**(1-real(r)/2)
1987  end if
1988  end if
1989  tcacc(r) = tcerr_aux(r)/norm(r)
1990  end do
1991 
1992  end if
1993 
1994  call propagateaccflag_cll(tcacc,rmax)
1995  call propagateerrflag_cll
1996 
1997  if (monitoring) then
1998  pointscntcten_cll = pointscntcten_cll + 1
1999 
2000  if(maxval(tcacc).gt.reqacc_cll) accpointscntcten_cll = accpointscntcten_cll + 1
2001 
2002  if(maxval(tcacc).gt.critacc_cll) then
2003  critpointscntcten_cll = critpointscntcten_cll + 1
2004  if ( critpointscntcten_cll.le.noutcritpointsmax_cll(3) ) then
2005  call critpointsout_cll('TCten_cll',0,maxval(tcacc),critpointscntcten_cll)
2006  if( critpointscntcten_cll.eq.noutcritpointsmax_cll(3)) then
2007  write(ncpout_cll,*) ' Further output of Critical Points for TCten_cll suppressed'
2008  write(ncpout_cll,*)
2009  endif
2010 #ifdef CritPoints2
2011  call critpointsout2_cll('TCten_cll',0,maxval(tcacc),critpointscntcten_cll)
2012  if( critpointscntcten_cll.eq.noutcritpointsmax_cll(3)) then
2013  write(ncpout2_cll,*) ' Further output of Critical Points for TCten_cll suppressed'
2014  write(ncpout2_cll,*)
2015  endif
2016 #endif
2017  end if
2018  end if
2019  end if
2020 

◆ cten_args_list_cll()

subroutine collier_tensors::cten_cll::cten_args_list_cll ( double complex, dimension(:), intent(out)  TC,
double complex, dimension(:), intent(out)  TCuv,
double complex, dimension(0:3), intent(in)  p1vec,
double complex, dimension(0:3), intent(in)  p2vec,
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  TCerr 
)

Definition at line 2033 of file collier_tensors.F90.

2033  integer, intent(in) :: rmax
2034  double complex, intent(in) :: p1vec(0:3), p2vec(0:3)
2035  double complex, intent(in) :: p10,p21,p20,m02,m12,m22
2036  double complex, intent(out) :: TC(:), TCuv(:)
2037  double precision, intent(out), optional :: TCerr(0:rmax)
2038  logical :: eflag
2039 
2040  if (3.gt.nmax_cll) then
2041  call seterrflag_cll(-10)
2042  call errout_cll('Cten_cll','Nmax_cll smaller 3',eflag,.true.)
2043  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
2044  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 3'
2045  call propagateerrflag_cll
2046  return
2047  end if
2048  if (rmax.gt.rmax_cll) then
2049  call seterrflag_cll(-10)
2050  call errout_cll('Cten_cll','argument rmax larger than rmax_cll',eflag,.true.)
2051  write(nerrout_cll,*) 'rmax =',rmax
2052  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
2053  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
2054  call propagateerrflag_cll
2055  return
2056  end if
2057 
2058  call cten_args_list_checked_cll(tc,tcuv,p1vec,p2vec,p10,p21,p20,m02,m12,m22,rmax,tcerr)
2059 

◆ cten_list_cll()

subroutine collier_tensors::cten_cll::cten_list_cll ( double complex, dimension(:), intent(out)  TC,
double complex, dimension(:), intent(out)  TCuv,
double complex, dimension(0:3,2), intent(in)  MomVec,
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  TCerr 
)

Definition at line 1669 of file collier_tensors.F90.

1669 
1670  integer, intent(in) :: rmax
1671  double complex, intent(in) :: MomVec(0:3,2), MomInv(3), masses2(0:2)
1672  double complex, intent(out) :: TC(:), TCuv(:)
1673  double precision, intent(out), optional :: TCerr(0:rmax)
1674  logical :: eflag
1675 
1676  if (3.gt.nmax_cll) then
1677  call seterrflag_cll(-10)
1678  call errout_cll('Cten_cll','Nmax_cll smaller 3',eflag,.true.)
1679  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
1680  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 3'
1681  call propagateerrflag_cll
1682  return
1683  end if
1684  if (rmax.gt.rmax_cll) then
1685  call seterrflag_cll(-10)
1686  call errout_cll('Cten_cll','argument rmax larger than rmax_cll',eflag,.true.)
1687  write(nerrout_cll,*) 'rmax =',rmax
1688  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
1689  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
1690  call propagateerrflag_cll
1691  return
1692  end if
1693 
1694  call cten_list_checked_cll(tc,tcuv,momvec,mominv,masses2,rmax,tcerr)
1695 

◆ cten_main_cll()

subroutine collier_tensors::cten_cll::cten_main_cll ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TC,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TCuv,
double complex, dimension(0:3,2), intent(in)  MomVec,
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  TCerr 
)

Definition at line 1494 of file collier_tensors.F90.

1494 
1495  integer, intent(in) :: rmax
1496  double complex, intent(in) :: MomVec(0:3,2), MomInv(3), masses2(0:2)
1497  double complex, intent(out) :: TC(0:rmax,0:rmax,0:rmax,0:rmax)
1498  double complex, intent(out) :: TCuv(0:rmax,0:rmax,0:rmax,0:rmax)
1499  double precision, intent(out), optional :: TCerr(0:rmax)
1500  double complex :: TC2(0:rmax,0:rmax,0:rmax,0:rmax), TCuv2(0:rmax,0:rmax,0:rmax,0:rmax)
1501  double complex :: CC(0:rmax/2,0:rmax,0:rmax), CCuv(0:rmax/2,0:rmax,0:rmax)
1502  double precision :: CCerr(0:rmax), TCerr_aux(0:rmax), TCerr_aux2(0:rmax), TCacc(0:rmax)
1503  double complex args(14)
1504  double precision :: TCdiff(0:rmax),norm(0:rmax),norm_coli,norm_dd
1505  integer :: r,n0,n1,n2,n3
1506  logical :: eflag
1507 
1508  if (3.gt.nmax_cll) then
1509  call seterrflag_cll(-10)
1510  call errout_cll('Cten_cll','Nmax_cll smaller 3',eflag,.true.)
1511  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
1512  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 3'
1513  call propagateerrflag_cll
1514  return
1515  end if
1516  if (rmax.gt.rmax_cll) then
1517  call seterrflag_cll(-10)
1518  call errout_cll('Cten_cll','argument rmax larger than rmax_cll',eflag,.true.)
1519  write(nerrout_cll,*) 'rmax =',rmax
1520  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
1521  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
1522  call propagateerrflag_cll
1523  return
1524  end if
1525 
1526  ! set ID of master call
1527  args(1:4) = momvec(0:,1)
1528  args(5:8) = momvec(0:,2)
1529  args(9:11) = mominv
1530  args(12:14) = masses2(0:)
1531  call setmasterfname_cll('Cten_cll')
1532  call setmastern_cll(3)
1533  call setmasterr_cll(rmax)
1534  call setmasterargs_cll(14,args)
1535 
1536  call settencache_cll(tenred_cll-1)
1537 
1538  if (mode_cll.eq.3) then
1539  ! calculate tensor with coefficients from COLI
1540  mode_cll = 1
1541  call c_main_cll(cc,ccuv,mominv(1),mominv(2),mominv(3), &
1542  masses2(0),masses2(1),masses2(2),rmax,cerr2=ccerr,id_in=0)
1543  call calctensorc(tc,tcuv,tcerr_aux,cc,ccuv,ccerr,momvec,rmax)
1544 
1545  ! calculate tensor with coefficients from DD
1546  mode_cll = 2
1547  call c_main_cll(cc,ccuv,mominv(1),mominv(2),mominv(3), &
1548  masses2(0),masses2(1),masses2(2),rmax,cerr2=ccerr,id_in=0)
1549  call calctensorc(tc2,tcuv2,tcerr_aux2,cc,ccuv,ccerr,momvec,rmax)
1550 
1551  ! comparison --> take better result
1552  mode_cll = 3
1553  do r=0,rmax
1554  norm_coli=0d0
1555  norm_dd=0d0
1556  do n0=0,r
1557  do n1=0,r-n0
1558  do n2=0,r-n0-n1
1559  n3=r-n0-n1-n2
1560  norm_coli = max(norm_coli,abs(tc(n0,n1,n2,n3)))
1561  norm_dd = max(norm_dd,abs(tc2(n0,n1,n2,n3)))
1562  end do
1563  end do
1564  end do
1565  if (norm_coli.eq.0d0) then
1566  norm_coli = max(maxval(abs(mominv(1:3))),maxval(abs(masses2(0:2))))
1567  if(norm_coli.ne.0d0) then
1568  norm_coli=1d0/norm_coli**(1-real(r)/2)
1569  else
1570  norm_coli=1d0/muir2_cll**(1-real(r)/2)
1571  end if
1572  end if
1573  if (norm_dd.eq.0d0) then
1574  norm_dd = max(maxval(abs(mominv(1:3))),maxval(abs(masses2(0:2))))
1575  if(norm_dd.ne.0d0) then
1576  norm_dd=1d0/norm_dd**(1-real(r)/2)
1577  else
1578  norm_dd=1d0/muir2_cll**(1-real(r)/2)
1579  end if
1580  end if
1581  norm(r) = min(norm_coli,norm_dd)
1582  end do
1583 
1584  call checktensors_cll(tc,tc2,momvec,mominv,masses2,norm,3,rmax,tcdiff)
1585 
1586  if (tcerr_aux(rmax).lt.tcerr_aux2(rmax)) then
1587  if (present(tcerr)) tcerr = max(tcerr_aux,tcdiff*norm)
1588  do r=0,rmax
1589  tcacc(r) = max(tcerr_aux(r)/norm(r),tcdiff(r))
1590  end do
1591  if (monitoring) pointscntcten_coli = pointscntcten_coli + 1
1592  else
1593  tc = tc2
1594  tcuv = tcuv2
1595  if (present(tcerr)) tcerr = max(tcerr_aux2,tcdiff*norm)
1596  do r=0,rmax
1597  tcacc(r) = max(tcerr_aux2(r)/norm(r),tcdiff(r))
1598  end do
1599  if (monitoring) pointscntcten_dd = pointscntcten_dd + 1
1600  end if
1601 
1602  else
1603  call c_main_cll(cc,ccuv,mominv(1),mominv(2),mominv(3), &
1604  masses2(0),masses2(1),masses2(2),rmax,cerr2=ccerr,id_in=0)
1605  call calctensorc(tc,tcuv,tcerr_aux,cc,ccuv,ccerr,momvec,rmax)
1606  if (present(tcerr)) tcerr = tcerr_aux
1607  norm=0d0
1608  do r=0,rmax
1609  do n0=0,r
1610  do n1=0,r-n0
1611  do n2=0,r-n0-n1
1612  n3=r-n0-n1-n2
1613  norm(r) = max(norm(r),abs(tc(n0,n1,n2,n3)))
1614  end do
1615  end do
1616  end do
1617  if (norm(r).eq.0d0) then
1618  norm(r) = max(maxval(abs(mominv(1:3))),maxval(abs(masses2(0:2))))
1619  if(norm(r).ne.0d0) then
1620  norm(r)=1d0/norm(r)**(1-real(r)/2)
1621  else
1622  norm(r)=1d0/muir2_cll**(1-real(r)/2)
1623  end if
1624  end if
1625  tcacc(r) = tcerr_aux(r)/norm(r)
1626  end do
1627 
1628  end if
1629 
1630  call propagateaccflag_cll(tcacc,rmax)
1631  call propagateerrflag_cll
1632 
1633  if (monitoring) then
1634  pointscntcten_cll = pointscntcten_cll + 1
1635 
1636  if(maxval(tcacc).gt.reqacc_cll) accpointscntcten_cll = accpointscntcten_cll + 1
1637 
1638  if(maxval(tcacc).gt.critacc_cll) then
1639  critpointscntcten_cll = critpointscntcten_cll + 1
1640  if ( critpointscntcten_cll.le.noutcritpointsmax_cll(3) ) then
1641  call critpointsout_cll('TCten_cll',0,maxval(tcacc),critpointscntcten_cll)
1642  if( critpointscntcten_cll.eq.noutcritpointsmax_cll(3)) then
1643  write(ncpout_cll,*) ' Further output of Critical Points for TCten_cll suppressed'
1644  write(ncpout_cll,*)
1645  endif
1646 #ifdef CritPoints2
1647  call critpointsout2_cll('TCten_cll',0,maxval(tcacc),critpointscntcten_cll)
1648  if( critpointscntcten_cll.eq.noutcritpointsmax_cll(3)) then
1649  write(ncpout2_cll,*) ' Further output of Critical Points for TCten_cll suppressed'
1650  write(ncpout2_cll,*)
1651  endif
1652 #endif
1653  end if
1654  end if
1655  end if
1656 

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