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::e_cll Interface Reference

Public Member Functions

subroutine e_main_cll (E, Euv, p10, p21, p32, p43, p40, p20, p31, p42, p30, p41, m02, m12, m22, m32, m42, rmax, Eerr, id_in, Eerr2)
 
subroutine e_arrays_cll (E, Euv, MomInv, masses2, rmax, Eerr, Eerr2)
 
subroutine e_list_cll (E, Euv, p10, p21, p32, p43, p40, p20, p31, p42, p30, p41, m02, m12, m22, m32, m42, rmax, Eerr, Eerr2)
 
subroutine e_arrays_list_cll (E, Euv, MomInv, masses2, rmax, Eerr, Eerr2)
 

Detailed Description

Definition at line 58 of file collier_coefs.F90.

Member Function/Subroutine Documentation

◆ e_arrays_cll()

subroutine collier_coefs::e_cll::e_arrays_cll ( double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  E,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  Euv,
double complex, dimension(10), intent(in)  MomInv,
double complex, dimension(0:4), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Eerr,
double precision, dimension(0:rmax), intent(out), optional  Eerr2 
)

Definition at line 2137 of file collier_coefs.F90.

2137 
2138  integer, intent(in) :: rmax
2139  double complex, intent(in) :: MomInv(10), masses2(0:4)
2140  double complex, intent(out) :: E(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
2141  double complex, intent(out) :: Euv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
2142  double precision, optional, intent(out) :: Eerr(0:rmax),Eerr2(0:rmax)
2143  double precision :: Eerraux(0:rmax),Eerr2aux(0:rmax)
2144 
2145  if (present(eerr)) then
2146  if (present(eerr2)) then
2147  call e_main_cll(e,euv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5), &
2148  mominv(6),mominv(7),mominv(8),mominv(9),mominv(10), &
2149  masses2(0),masses2(1),masses2(2),masses2(3),masses2(4),rmax,eerr,eerr2=eerr2)
2150  else
2151  call e_main_cll(e,euv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5), &
2152  mominv(6),mominv(7),mominv(8),mominv(9),mominv(10), &
2153  masses2(0),masses2(1),masses2(2),masses2(3),masses2(4),rmax,eerr)
2154  end if
2155  else
2156  if (present(eerr2)) then
2157  call e_main_cll(e,euv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5), &
2158  mominv(6),mominv(7),mominv(8),mominv(9),mominv(10), &
2159  masses2(0),masses2(1),masses2(2),masses2(3),masses2(4),rmax,eerraux,eerr2=eerr2)
2160  else
2161  call e_main_cll(e,euv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5), &
2162  mominv(6),mominv(7),mominv(8),mominv(9),mominv(10), &
2163  masses2(0),masses2(1),masses2(2),masses2(3),masses2(4),rmax,eerraux)
2164  end if
2165  end if
2166 

◆ e_arrays_list_cll()

subroutine collier_coefs::e_cll::e_arrays_list_cll ( double complex, dimension(:), intent(out)  E,
double complex, dimension(:), intent(out)  Euv,
double complex, dimension(10), intent(in)  MomInv,
double complex, dimension(0:4), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Eerr,
double precision, dimension(0:rmax), intent(out), optional  Eerr2 
)

Definition at line 2276 of file collier_coefs.F90.

2276 
2277  integer, intent(in) :: rmax
2278  double complex, intent(in) :: MomInv(10), masses2(0:4)
2279  double complex, intent(out) :: E(:),Euv(:)
2280  double precision, optional, intent(out) :: Eerr(0:rmax),Eerr2(0:rmax)
2281  logical :: eflag
2282 
2283  if (5.gt.nmax_cll) then
2284  call seterrflag_cll(-10)
2285  call errout_cll('E_cll','Nmax_cll smaller 5',eflag,.true.)
2286  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
2287  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 5'
2288  call propagateerrflag_cll
2289  return
2290  end if
2291  if (rmax.gt.rmax_cll) then
2292  call seterrflag_cll(-10)
2293  call errout_cll('E_cll','argument rmax larger than rmax_cll',eflag,.true.)
2294  write(nerrout_cll,*) 'rmax =',rmax
2295  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
2296  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
2297  call propagateerrflag_cll
2298  return
2299  end if
2300 
2301  call e_arrays_list_checked_cll(e,euv,mominv,masses2,rmax,eerr,eerr2)
2302 

◆ e_list_cll()

subroutine collier_coefs::e_cll::e_list_cll ( double complex, dimension(:), intent(out)  E,
double complex, dimension(:), intent(out)  Euv,
double complex, intent(in)  p10,
double complex, intent(in)  p21,
double complex, intent(in)  p32,
double complex, intent(in)  p43,
double complex, intent(in)  p40,
double complex, intent(in)  p20,
double complex, intent(in)  p31,
double complex, intent(in)  p42,
double complex, intent(in)  p30,
double complex, intent(in)  p41,
double complex, intent(in)  m02,
double complex, intent(in)  m12,
double complex, intent(in)  m22,
double complex, intent(in)  m32,
double complex, intent(in)  m42,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Eerr,
double precision, dimension(0:rmax), intent(out), optional  Eerr2 
)

Definition at line 2182 of file collier_coefs.F90.

2182 
2183  integer, intent(in) :: rmax
2184  double complex, intent(in) :: p10,p21,p32,p43,p40,p20,p31,p42,p30,p41
2185  double complex, intent(in) :: m02,m12,m22,m32,m42
2186  double complex, intent(out) :: E(:),Euv(:)
2187  double precision, optional, intent(out) :: Eerr(0:rmax),Eerr2(0:rmax)
2188  logical :: eflag
2189 
2190  if (5.gt.nmax_cll) then
2191  call seterrflag_cll(-10)
2192  call errout_cll('E_cll','Nmax_cll smaller 5',eflag,.true.)
2193  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
2194  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 5'
2195  call propagateerrflag_cll
2196  return
2197  end if
2198  if (rmax.gt.rmax_cll) then
2199  call seterrflag_cll(-10)
2200  call errout_cll('E_cll','argument rmax larger than rmax_cll',eflag,.true.)
2201  write(nerrout_cll,*) 'rmax =',rmax
2202  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
2203  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
2204  call propagateerrflag_cll
2205  return
2206  end if
2207 
2208  call e_list_checked_cll(e,euv,p10,p21,p32,p43,p40,p20,p31,p42,p30,p41, &
2209  m02,m12,m22,m32,m42,rmax,eerr,eerr2)
2210 

◆ e_main_cll()

subroutine collier_coefs::e_cll::e_main_cll ( double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  E,
double complex, dimension(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  Euv,
double complex, intent(in)  p10,
double complex, intent(in)  p21,
double complex, intent(in)  p32,
double complex, intent(in)  p43,
double complex, intent(in)  p40,
double complex, intent(in)  p20,
double complex, intent(in)  p31,
double complex, intent(in)  p42,
double complex, intent(in)  p30,
double complex, intent(in)  p41,
double complex, intent(in)  m02,
double complex, intent(in)  m12,
double complex, intent(in)  m22,
double complex, intent(in)  m32,
double complex, intent(in)  m42,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  Eerr,
integer, intent(in), optional  id_in,
double precision, dimension(0:rmax), intent(out), optional  Eerr2 
)

Definition at line 1780 of file collier_coefs.F90.

1780 
1781  integer, intent(in) :: rmax
1782  double complex, intent(in) :: p10,p21,p32,p43,p40,p20,p31,p42,p30,p41
1783  double complex, intent(in) :: m02,m12,m22,m32,m42
1784  double complex, intent(out) :: E(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
1785  double complex, intent(out) :: Euv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
1786  double precision, optional, intent(out) :: Eerr(0:rmax),Eerr2(0:rmax)
1787  double precision :: q10,q21,q32,q43,q40,q20,q31,q42,q30,q41
1788  double complex :: mm02,mm12,mm22,mm32,mm42
1789  double precision :: Eerraux(0:rmax),Eerr2aux(0:rmax),Ediff(0:rmax)
1790  integer, optional, intent(in) :: id_in
1791  double complex :: E2uv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
1792  double complex :: E2(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
1793  double complex :: Edd(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
1794  double complex :: elimcminf2
1795  double complex :: args(15)
1796  integer :: n0,rank,errflag,id
1797  double precision :: accrelDD(0:rmax_DD),accabsDD(0:rmax_DD),Eacc(0:rmax),norm,norm_coli,norm_dd,Eacc2(0:rmax)
1798  double precision :: accrel2DD(0:rmax_DD),accabs2DD(0:rmax_DD)
1799  integer :: accflagDD,errflagDD,NDD,rankDD
1800  logical :: mflag,eflag
1801  integer :: r,n1,n2,n3,n4
1802 
1803  if (5.gt.nmax_cll) then
1804  call seterrflag_cll(-10)
1805  call errout_cll('E_cll','Nmax_cll smaller 5',eflag,.true.)
1806  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
1807  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 5'
1808  call propagateerrflag_cll
1809  return
1810  end if
1811  if (rmax.gt.rmax_cll) then
1812  call seterrflag_cll(-10)
1813  call errout_cll('E_cll','argument rmax larger than rmax_cll',eflag,.true.)
1814  write(nerrout_cll,*) 'rmax =',rmax
1815  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
1816  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
1817  call propagateerrflag_cll
1818  return
1819  end if
1820 
1821  mflag=.true.
1822  if (present(id_in)) then
1823  mflag=.false.
1824  id = id_in
1825  else
1826  id = 0
1827  end if
1828 
1829  if (mflag) then
1830  ! set ID of master call
1831  args(1) = p10
1832  args(2) = p21
1833  args(3) = p32
1834  args(4) = p43
1835  args(5) = p40
1836  args(6) = p20
1837  args(7) = p31
1838  args(8) = p42
1839  args(9) = p30
1840  args(10) = p41
1841  args(11) = m02
1842  args(12) = m12
1843  args(13) = m22
1844  args(14) = m32
1845  args(15) = m42
1846  call setmasterfname_cll('E_cll')
1847  call setmastern_cll(5)
1848  call setmasterr_cll(rmax)
1849  call setmasterargs_cll(15,args)
1850 
1851  call settencache_cll(never_tenred_cll)
1852  end if
1853 
1854 
1855  select case (mode_cll)
1856 
1857  case (1)
1858  ! calculate loop integral using
1859  ! COLI implementation by AD/LH
1860 
1861  call calce(e,euv,p10,p21,p32,p43,p40,p20,p31,p42,p30,p41, &
1862  m02,m12,m22,m32,m42,rmax,id,eerraux,eerr2aux)
1863 
1864  norm = abs(e(0,0,0,0,0))
1865  do r=1,rmax
1866  do n1=0,r
1867  do n2=0,r-n1
1868  do n3=0,r-n1-n2
1869  n4=r-n1-n2-n3
1870  norm = max(norm,abs(e(0,n1,n2,n3,n4)))
1871  end do
1872  end do
1873  end do
1874  end do
1875  if (norm.eq.0d0) then
1876  norm = max(abs(p10),abs(p21),abs(p32),abs(p43),abs(p40), &
1877  abs(p20),abs(p31),abs(p42),abs(p30),abs(p41), &
1878  abs(m02),abs(m12),abs(m22),abs(m32),abs(m42))
1879  if(norm.ne.0d0) then
1880  norm=1d0/norm**3
1881  else
1882  norm=1d0/muir2_cll**3
1883  end if
1884  end if
1885  if (norm.ne.0d0) then
1886  eacc = eerraux/norm
1887  eacc2 = eerr2aux/norm
1888  else
1889  eacc = 0d0
1890  eacc2 = 0d0
1891  end if
1892 
1893  if (present(eerr)) eerr = eerraux
1894  if (present(eerr2)) eerr2 = eerr2aux
1895 
1896  if (mflag) call propagateaccflag_cll(eacc,rmax)
1897 
1898 
1899  case (2)
1900  ! calculate loop integral using
1901  ! DD implementation by SD
1902 
1903  if (rmax.gt.5) then
1904  call seterrflag_cll(-10)
1905  call errout_cll('E_cll','rank higher than maximum rank implemented in DD library',eflag)
1906  if(eflag) then
1907  write(nerrout_cll,*) 'E_cll: 5-point function of rank>5 not implemented in DD library'
1908  end if
1909  end if
1910 
1911 
1912  ! replace small masses by DD-identifiers
1913  q10 = dreal(getminf2dd_cll(p10))
1914  q21 = dreal(getminf2dd_cll(p21))
1915  q32 = dreal(getminf2dd_cll(p32))
1916  q43 = dreal(getminf2dd_cll(p43))
1917  q40 = dreal(getminf2dd_cll(p40))
1918  q20 = dreal(getminf2dd_cll(p20))
1919  q31 = dreal(getminf2dd_cll(p31))
1920  q42 = dreal(getminf2dd_cll(p42))
1921  q30 = dreal(getminf2dd_cll(p30))
1922  q41 = dreal(getminf2dd_cll(p41))
1923  mm02 = getminf2dd_cll(m02)
1924  mm12 = getminf2dd_cll(m12)
1925  mm22 = getminf2dd_cll(m22)
1926  mm32 = getminf2dd_cll(m32)
1927  mm42 = getminf2dd_cll(m42)
1928 
1929  rank = rmax
1930  call e_dd(edd,q10,q21,q32,q43,q40,q20,q31,q42,q30,q41, &
1931  mm02,mm12,mm22,mm32,mm42,rank,id)
1932  e(0:rank/2,0:rank,0:rank,0:rank,0:rank) = edd(0:rank/2,0:rank,0:rank,0:rank,0:rank)
1933  euv = 0d0
1934 
1935  call ddgetacc(accreldd,accabsdd,accrel2dd,accabs2dd,ndd,rankdd,accflagdd,errflagdd,id)
1936  if (present(eerr)) eerr(0:rmax) = accabsdd(0:rmax)
1937  if (present(eerr2)) eerr2(0:rmax) = accabs2dd(0:rmax)
1938 
1939  norm = abs(e(0,0,0,0,0))
1940  do r=1,rmax
1941  do n1=0,r
1942  do n2=0,r-n1
1943  do n3=0,r-n1-n2
1944  n4=r-n1-n2-n3
1945  norm = max(norm,abs(e(0,n1,n2,n3,n4)))
1946  end do
1947  end do
1948  end do
1949  end do
1950  if (norm.eq.0d0) then
1951  norm = max(abs(p10),abs(p21),abs(p32),abs(p43),abs(p40), &
1952  abs(p20),abs(p31),abs(p42),abs(p30),abs(p41), &
1953  abs(m02),abs(m12),abs(m22),abs(m32),abs(m42))
1954  if(norm.ne.0d0) then
1955  norm=1d0/norm**3
1956  else
1957  norm=1d0/muir2_cll**3
1958  end if
1959  end if
1960  if (norm.ne.0d0) then
1961  eacc = accabsdd(0:rmax)/norm
1962  eacc2 = accabs2dd(0:rmax)/norm
1963  else
1964  eacc = 0d0
1965  eacc2 = 0d0
1966  end if
1967  if (mflag) call propagateaccflag_cll(eacc,rmax)
1968 
1969 
1970  case (3)
1971  ! cross-check mode
1972  ! compare results for loop integral
1973 
1974 
1975  ! from COLI implementation by AD/LH and
1976  ! from DD implementation by SD
1977 
1978  ! calculate loop integral
1979  call calce(e,euv,p10,p21,p32,p43,p40,p20,p31,p42,p30,p41, &
1980  m02,m12,m22,m32,m42,rmax,id,eerraux,eerr2aux)
1981 
1982 
1983  ! calculate loop integral
1984  if (rmax.gt.5) then
1985  call seterrflag_cll(-10)
1986  call errout_cll('E_cll','rank higher than maximum rank implemented in DD library',eflag)
1987  if(eflag) then
1988  write(nerrout_cll,*) 'E_cll: 5-point function of rank>5 not implemented in DD library'
1989  end if
1990  end if
1991 
1992 
1993  ! replace small masses by DD-identifiers
1994  q10 = dreal(getminf2dd_cll(p10))
1995  q21 = dreal(getminf2dd_cll(p21))
1996  q32 = dreal(getminf2dd_cll(p32))
1997  q43 = dreal(getminf2dd_cll(p43))
1998  q40 = dreal(getminf2dd_cll(p40))
1999  q20 = dreal(getminf2dd_cll(p20))
2000  q31 = dreal(getminf2dd_cll(p31))
2001  q42 = dreal(getminf2dd_cll(p42))
2002  q30 = dreal(getminf2dd_cll(p30))
2003  q41 = dreal(getminf2dd_cll(p41))
2004  mm02 = getminf2dd_cll(m02)
2005  mm12 = getminf2dd_cll(m12)
2006  mm22 = getminf2dd_cll(m22)
2007  mm32 = getminf2dd_cll(m32)
2008  mm42 = getminf2dd_cll(m42)
2009 
2010  rank = rmax
2011  call e_dd(edd,q10,q21,q32,q43,q40,q20,q31,q42,q30,q41, &
2012  mm02,mm12,mm22,mm32,mm42,rank,id)
2013  e2(0:rank/2,0:rank,0:rank,0:rank,0:rank) = edd(0:rank/2,0:rank,0:rank,0:rank,0:rank)
2014  e2uv = 0d0
2015 
2016  call ddgetacc(accreldd,accabsdd,accrel2dd,accabs2dd,ndd,rankdd,accflagdd,errflagdd,id)
2017 
2018  norm_coli = abs(e(0,0,0,0,0))
2019  norm_dd = abs(e2(0,0,0,0,0))
2020  do r=1,rmax
2021  do n1=0,r
2022  do n2=0,r-n1
2023  do n3=0,r-n1-n2
2024  n4=r-n1-n2-n3
2025  norm_coli = max(norm_coli,abs(e(0,n1,n2,n3,n4)))
2026  norm_dd = max(norm_dd,abs(e2(0,n1,n2,n3,n4)))
2027  end do
2028  end do
2029  end do
2030  end do
2031  if (norm_coli.eq.0d0) then
2032  norm_coli = max(abs(p10),abs(p21),abs(p32),abs(p43),abs(p40), &
2033  abs(p20),abs(p31),abs(p42),abs(p30),abs(p41), &
2034  abs(m02),abs(m12),abs(m22),abs(m32),abs(m42))
2035  if(norm_coli.ne.0d0) then
2036  norm_coli=1d0/norm_coli**3
2037  else
2038  norm_coli=1d0/muir2_cll**3
2039  end if
2040  end if
2041  if (norm_dd.eq.0d0) then
2042  norm_dd = max(abs(p10),abs(p21),abs(p32),abs(p43),abs(p40), &
2043  abs(p20),abs(p31),abs(p42),abs(p30),abs(p41), &
2044  abs(m02),abs(m12),abs(m22),abs(m32),abs(m42))
2045  if(norm_dd.ne.0d0) then
2046  norm_dd=1d0/norm_dd**3
2047  else
2048  norm_dd=1d0/muir2_cll**3
2049  end if
2050  end if
2051  norm=min(norm_coli,norm_dd)
2052 
2053  ! cross-check
2054  call checkcoefse_cll(e,e2,p10,p21,p32,p43,p40,p20,p31,p42,p30,p41, &
2055  m02,m12,m22,m32,m42,rmax,norm,ediff)
2056 
2057 
2058  if (eerraux(rmax).lt.accabsdd(rmax)) then
2059  if (present(eerr)) eerr = max(eerraux,ediff)
2060  if (present(eerr2)) eerr2 = eerr2aux
2061  if (norm.ne.0d0) then
2062  eacc = max(eerraux/norm_coli,ediff/norm)
2063  eacc2 = eerr2aux/norm_coli
2064  else
2065  eacc = ediff
2066  eacc2 = 0d0
2067  end if
2068  if (monitoring) pointscnte_coli = pointscnte_coli + 1
2069  else
2070  e = e2
2071  euv = e2uv
2072  if (present(eerr)) eerr = max(accabsdd(0:rmax),ediff)
2073  if (present(eerr2)) eerr2 = accabs2dd(0:rmax)
2074  if (norm.ne.0d0) then
2075  eacc = max(accabsdd(0:rmax)/norm_dd,ediff/norm)
2076  eacc2 = accabs2dd(0:rmax)/norm_dd
2077  else
2078  eacc = ediff
2079  eacc2 = 0d0
2080  end if
2081  if (monitoring) pointscnte_dd = pointscnte_dd + 1
2082  end if
2083 
2084  if (mflag) call propagateaccflag_cll(eacc,rmax)
2085 
2086  end select
2087 
2088  if (mflag) call propagateerrflag_cll
2089 
2090  if (monitoring) then
2091  pointscnte_cll = pointscnte_cll + 1
2092 
2093  if(maxval(eacc).gt.reqacc_cll) accpointscnte_cll = accpointscnte_cll + 1
2094 
2095  if(maxval(eacc).gt.critacc_cll) then
2096  critpointscnte_cll = critpointscnte_cll + 1
2097  if ( critpointscnte_cll.le.noutcritpointsmax_cll(5) ) then
2098  call critpointsout_cll('E_cll',0,maxval(eacc), critpointscnte_cll)
2099  if( critpointscnte_cll.eq.noutcritpointsmax_cll(5)) then
2100  write(ncpout_cll,*) ' Further output of Critical Points for E_cll suppressed '
2101  write(ncpout_cll,*)
2102  endif
2103  end if
2104  end if
2105 
2106 #ifdef CritPoints2
2107  if(maxval(eacc2).gt.reqacc_cll) accpointscnte2_cll = accpointscnte2_cll + 1
2108 
2109  if(maxval(eacc2).gt.critacc_cll) then
2110  critpointscnte2_cll = critpointscnte2_cll + 1
2111  if ( critpointscnte2_cll.le.noutcritpointsmax_cll(5) ) then
2112  call critpointsout2_cll('E_cll',0,maxval(eacc2), critpointscnte2_cll)
2113  if( critpointscnte2_cll.eq.noutcritpointsmax_cll(5)) then
2114  write(ncpout2_cll,*) ' Further output of Critical Points for E_cll suppressed '
2115  write(ncpout2_cll,*)
2116  endif
2117  end if
2118  end if
2119 #endif
2120 
2121  end if
2122 
2123 

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