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_tensors::dten_cll Interface Reference

Public Member Functions

subroutine dten_main_cll (TD, TDuv, MomVec, MomInv, masses2, rmax, TDerr)
 
subroutine dten_list_cll (TD, TDuv, MomVec, MomInv, masses2, rmax, TDerr)
 
subroutine dten_args_cll (TD, TDuv, p1vec, p2vec, p3vec, p10, p21, p32, p30, p20, p31, m02, m12, m22, m32, rmax, TDerr)
 
subroutine dten_args_list_cll (TD, TDuv, p1vec, p2vec, p3vec, p10, p21, p32, p30, p20, p31, m02, m12, m22, m32, rmax, TDerr)
 

Detailed Description

Definition at line 62 of file collier_tensors.F90.

Member Function/Subroutine Documentation

◆ dten_args_cll()

subroutine collier_tensors::dten_cll::dten_args_cll ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TD,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TDuv,
double complex, dimension(0:3), intent(in)  p1vec,
double complex, dimension(0:3), intent(in)  p2vec,
double complex, dimension(0:3), intent(in)  p3vec,
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  TDerr 
)

Definition at line 2580 of file collier_tensors.F90.

2580 
2581  integer, intent(in) :: rmax
2582  double complex, intent(in) :: p1vec(0:3), p2vec(0:3), p3vec(0:3)
2583  double complex, intent(in) :: p10,p21,p32,p30,p20,p31,m02,m12,m22,m32
2584  double complex, intent(out) :: TD(0:rmax,0:rmax,0:rmax,0:rmax)
2585  double complex, intent(out) :: TDuv(0:rmax,0:rmax,0:rmax,0:rmax)
2586  double precision, intent(out), optional :: TDerr(0:rmax)
2587  double complex TD2(0:rmax,0:rmax,0:rmax,0:rmax), TDuv2(0:rmax,0:rmax,0:rmax,0:rmax)
2588  double complex :: MomVec(0:3,3), MomInv(6), masses2(0:3)
2589  double complex :: CD(0:rmax/2,0:rmax,0:rmax,0:rmax)
2590  double complex :: CDuv(0:rmax/2,0:rmax,0:rmax,0:rmax)
2591  double precision :: CDerr(0:rmax), TDerr_aux(0:rmax), TDerr_aux2(0:rmax)
2592  double complex :: args(22)
2593  double precision :: TDdiff(0:rmax),norm(0:rmax),norm_coli,norm_dd,TDacc(0:rmax)
2594  integer :: r,n0,n1,n2,n3
2595  logical :: eflag
2596 
2597  if (4.gt.nmax_cll) then
2598  call seterrflag_cll(-10)
2599  call errout_cll('Dten_cll','Nmax_cll smaller 4',eflag,.true.)
2600  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
2601  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 4'
2602  call propagateerrflag_cll
2603  return
2604  end if
2605  if (rmax.gt.rmax_cll) then
2606  call seterrflag_cll(-10)
2607  call errout_cll('Dten_cll','argument rmax larger than rmax_cll',eflag,.true.)
2608  write(nerrout_cll,*) 'rmax =',rmax
2609  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
2610  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
2611  call propagateerrflag_cll
2612  return
2613  end if
2614 
2615  momvec(0:,1) = p1vec
2616  momvec(0:,2) = p2vec
2617  momvec(0:,3) = p3vec
2618  mominv(1) = p10
2619  mominv(2) = p21
2620  mominv(3) = p32
2621  mominv(4) = p30
2622  mominv(5) = p20
2623  mominv(6) = p31
2624  masses2(0) = m02
2625  masses2(1) = m12
2626  masses2(2) = m22
2627  masses2(3) = m32
2628 
2629  ! set ID of master call
2630  args(1:4) = momvec(0:,1)
2631  args(5:8) = momvec(0:,2)
2632  args(9:12) = momvec(0:,3)
2633  args(13:18) = mominv
2634  args(19:22) = masses2(0:)
2635  call setmasterfname_cll('Dten_cll')
2636  call setmastern_cll(4)
2637  call setmasterr_cll(rmax)
2638  call setmasterargs_cll(22,args)
2639 
2640  call settencache_cll(tenred_cll-1)
2641 
2642 
2643  if (mode_cll.eq.3) then
2644  ! calculate tensor with coefficients from COLI
2645  mode_cll = 1
2646  call d_main_cll(cd,cduv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
2647  masses2(0),masses2(1),masses2(2),masses2(3),rmax,derr2=cderr,id_in=0)
2648  call calctensord(td,tduv,tderr_aux,cd,cduv,cderr,momvec,rmax)
2649 
2650  ! calculate tensor with coefficients from DD
2651  mode_cll = 2
2652  call d_main_cll(cd,cduv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
2653  masses2(0),masses2(1),masses2(2),masses2(3),rmax,derr2=cderr,id_in=0)
2654  call calctensord(td2,tduv2,tderr_aux2,cd,cduv,cderr,momvec,rmax)
2655 
2656  ! comparison --> take better result
2657  mode_cll = 3
2658  do r=0,rmax
2659  norm_coli=0d0
2660  norm_dd=0d0
2661  do n0=0,r
2662  do n1=0,r-n0
2663  do n2=0,r-n0-n1
2664  n3=r-n0-n1-n2
2665  norm_coli = max(norm_coli,abs(td(n0,n1,n2,n3)))
2666  norm_dd = max(norm_dd,abs(td2(n0,n1,n2,n3)))
2667  end do
2668  end do
2669  end do
2670  if (norm_coli.eq.0d0) then
2671  norm_coli = max(maxval(abs(mominv(1:6))),maxval(abs(masses2(0:3))))
2672  if(norm_coli.ne.0d0) then
2673  norm_coli=1d0/norm_coli**(2-real(r)/2)
2674  else
2675  norm_coli=1d0/muir2_cll**(2-real(r)/2)
2676  end if
2677  end if
2678  if (norm_dd.eq.0d0) then
2679  norm_dd = max(maxval(abs(mominv(1:6))),maxval(abs(masses2(0:3))))
2680  if(norm_dd.ne.0d0) then
2681  norm_dd=1d0/norm_dd**(2-real(r)/2)
2682  else
2683  norm_dd=1d0/muir2_cll**(2-real(r)/2)
2684  end if
2685  end if
2686  norm(r) = min(norm_coli,norm_dd)
2687  end do
2688 
2689  call checktensors_cll(td,td2,momvec,mominv,masses2,norm,4,rmax,tddiff)
2690 
2691  if (tderr_aux(rmax).lt.tderr_aux2(rmax)) then
2692  if (present(tderr)) tderr = max(tderr_aux,tddiff*norm)
2693  do r=0,rmax
2694  tdacc(r) = max(tderr_aux(r)/norm(r),tddiff(r))
2695  end do
2696  if (monitoring) pointscntdten_coli = pointscntdten_coli + 1
2697  else
2698  td = td2
2699  tduv = tduv2
2700  if (present(tderr)) tderr = max(tderr_aux2,tddiff*norm)
2701  do r=0,rmax
2702  tdacc(r) = max(tderr_aux2(r)/norm(r),tddiff(r))
2703  end do
2704  if (monitoring) pointscntdten_dd = pointscntdten_dd + 1
2705  end if
2706 
2707  else
2708  call d_main_cll(cd,cduv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
2709  masses2(0),masses2(1),masses2(2),masses2(3),rmax,derr2=cderr,id_in=0)
2710 
2711  call calctensord(td,tduv,tderr_aux,cd,cduv,cderr,momvec,rmax)
2712  if (present(tderr)) tderr = tderr_aux
2713  norm=0d0
2714  do r=0,rmax
2715  do n0=0,r
2716  do n1=0,r-n0
2717  do n2=0,r-n0-n1
2718  n3=r-n0-n1-n2
2719  norm(r) = max(norm(r),abs(td(n0,n1,n2,n3)))
2720  end do
2721  end do
2722  end do
2723  if (norm(r).eq.0d0) then
2724  norm(r) = max(maxval(abs(mominv(1:6))),maxval(abs(masses2(0:3))))
2725  if(norm(r).ne.0d0) then
2726  norm(r)=1d0/norm(r)**(2-real(r)/2)
2727  else
2728  norm(r)=1d0/muir2_cll**(2-real(r)/2)
2729  end if
2730  end if
2731  tdacc(r) = tderr_aux(r)/norm(r)
2732  end do
2733 
2734  end if
2735 
2736  call propagateaccflag_cll(tdacc,rmax)
2737  call propagateerrflag_cll
2738 
2739  if (monitoring) then
2740  pointscntdten_cll = pointscntdten_cll + 1
2741 
2742  if(maxval(tdacc).gt.reqacc_cll) accpointscntdten_cll = accpointscntdten_cll + 1
2743 
2744  if(maxval(tdacc).gt.critacc_cll) then
2745  critpointscntdten_cll = critpointscntdten_cll + 1
2746  if ( critpointscntdten_cll.le.noutcritpointsmax_cll(4) ) then
2747  call critpointsout_cll('TDten_cll',0,maxval(tdacc),critpointscntdten_cll)
2748  if( critpointscntdten_cll.eq.noutcritpointsmax_cll(4)) then
2749  write(ncpout_cll,*) ' Further output of Critical Points for TDten_cll suppressed'
2750  write(ncpout_cll,*)
2751  endif
2752 #ifdef CritPoints2
2753  call critpointsout2_cll('TDten_cll',0,maxval(tdacc),critpointscntdten_cll)
2754  if( critpointscntdten_cll.eq.noutcritpointsmax_cll(4)) then
2755  write(ncpout2_cll,*) ' Further output of Critical Points for TDten_cll suppressed'
2756  write(ncpout2_cll,*)
2757  endif
2758 #endif
2759  end if
2760  end if
2761  end if
2762 

◆ dten_args_list_cll()

subroutine collier_tensors::dten_cll::dten_args_list_cll ( double complex, dimension(:), intent(out)  TD,
double complex, dimension(:), intent(out)  TDuv,
double complex, dimension(0:3), intent(in)  p1vec,
double complex, dimension(0:3), intent(in)  p2vec,
double complex, dimension(0:3), intent(in)  p3vec,
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  TDerr 
)

Definition at line 2777 of file collier_tensors.F90.

2777  integer, intent(in) :: rmax
2778  double complex, intent(in) :: p1vec(0:3), p2vec(0:3), p3vec(0:3)
2779  double complex, intent(in) :: p10,p21,p32,p30,p20,p31,m02,m12,m22,m32
2780  double complex, intent(out) :: TD(:), TDuv(:)
2781  double precision, intent(out), optional :: TDerr(0:rmax)
2782  logical :: eflag
2783 
2784  if (4.gt.nmax_cll) then
2785  call seterrflag_cll(-10)
2786  call errout_cll('Dten_cll','Nmax_cll smaller 4',eflag,.true.)
2787  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
2788  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 4'
2789  call propagateerrflag_cll
2790  return
2791  end if
2792  if (rmax.gt.rmax_cll) then
2793  call seterrflag_cll(-10)
2794  call errout_cll('Dten_cll','argument rmax larger than rmax_cll',eflag,.true.)
2795  write(nerrout_cll,*) 'rmax =',rmax
2796  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
2797  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
2798  call propagateerrflag_cll
2799  return
2800  end if
2801 
2802  call dten_args_list_checked_cll(td,tduv,p1vec,p2vec,p3vec,p10,p21,p32,p30,p20,p31, &
2803  m02,m12,m22,m32,rmax,tderr)
2804 

◆ dten_list_cll()

subroutine collier_tensors::dten_cll::dten_list_cll ( double complex, dimension(:), intent(out)  TD,
double complex, dimension(:), intent(out)  TDuv,
double complex, dimension(0:3,3), intent(in)  MomVec,
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  TDerr 
)

Definition at line 2400 of file collier_tensors.F90.

2400 
2401  integer, intent(in) :: rmax
2402  double complex, intent(in) :: MomVec(0:3,3), MomInv(6), masses2(0:3)
2403  double complex, intent(out) :: TD(:), TDuv(:)
2404  double precision, intent(out), optional :: TDerr(0:rmax)
2405  logical :: eflag
2406 
2407  if (4.gt.nmax_cll) then
2408  call seterrflag_cll(-10)
2409  call errout_cll('Dten_cll','Nmax_cll smaller 4',eflag,.true.)
2410  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
2411  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 4'
2412  call propagateerrflag_cll
2413  return
2414  end if
2415  if (rmax.gt.rmax_cll) then
2416  call seterrflag_cll(-10)
2417  call errout_cll('Dten_cll','argument rmax larger than rmax_cll',eflag,.true.)
2418  write(nerrout_cll,*) 'rmax =',rmax
2419  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
2420  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
2421  call propagateerrflag_cll
2422  return
2423  end if
2424 
2425  call dten_list_checked_cll(td,tduv,momvec,mominv,masses2,rmax,tderr)
2426 

◆ dten_main_cll()

subroutine collier_tensors::dten_cll::dten_main_cll ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TD,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TDuv,
double complex, dimension(0:3,3), intent(in)  MomVec,
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  TDerr 
)

Definition at line 2222 of file collier_tensors.F90.

2222 
2223  integer, intent(in) :: rmax
2224  double complex, intent(in) :: MomVec(0:3,3), MomInv(6), masses2(0:3)
2225  double complex, intent(out) :: TD(0:rmax,0:rmax,0:rmax,0:rmax)
2226  double complex, intent(out) :: TDuv(0:rmax,0:rmax,0:rmax,0:rmax)
2227  double precision, intent(out), optional :: TDerr(0:rmax)
2228  double complex :: CD(0:rmax/2,0:rmax,0:rmax,0:rmax)
2229  double complex :: TD2(0:rmax,0:rmax,0:rmax,0:rmax), TDuv2(0:rmax,0:rmax,0:rmax,0:rmax)
2230  double complex :: CDuv(0:rmax/2,0:rmax,0:rmax,0:rmax)
2231  double precision :: CDerr(0:rmax), TDerr_aux(0:rmax), TDerr_aux2(0:rmax)
2232  double complex :: args(22)
2233  double precision :: TDdiff(0:rmax),norm(0:rmax),norm_coli,norm_dd,TDacc(0:rmax)
2234  integer :: r,n0,n1,n2,n3
2235  logical :: eflag
2236 
2237  if (4.gt.nmax_cll) then
2238  call seterrflag_cll(-10)
2239  call errout_cll('Dten_cll','Nmax_cll smaller 4',eflag,.true.)
2240  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
2241  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 4'
2242  call propagateerrflag_cll
2243  return
2244  end if
2245  if (rmax.gt.rmax_cll) then
2246  call seterrflag_cll(-10)
2247  call errout_cll('Dten_cll','argument rmax larger than rmax_cll',eflag,.true.)
2248  write(nerrout_cll,*) 'rmax =',rmax
2249  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
2250  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
2251  call propagateerrflag_cll
2252  return
2253  end if
2254 
2255  ! set ID of master call
2256  args(1:4) = momvec(0:,1)
2257  args(5:8) = momvec(0:,2)
2258  args(9:12) = momvec(0:,3)
2259  args(13:18) = mominv
2260  args(19:22) = masses2(0:)
2261  call setmasterfname_cll('Dten_cll')
2262  call setmastern_cll(4)
2263  call setmasterr_cll(rmax)
2264  call setmasterargs_cll(22,args)
2265 
2266  call settencache_cll(tenred_cll-1)
2267 
2268  if (mode_cll.eq.3) then
2269  ! calculate tensor with coefficients from COLI
2270  mode_cll = 1
2271  call d_main_cll(cd,cduv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
2272  masses2(0),masses2(1),masses2(2),masses2(3),rmax,derr2=cderr,id_in=0)
2273  call calctensord(td,tduv,tderr_aux,cd,cduv,cderr,momvec,rmax)
2274 
2275  ! calculate tensor with coefficients from DD
2276  mode_cll = 2
2277  call d_main_cll(cd,cduv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
2278  masses2(0),masses2(1),masses2(2),masses2(3),rmax,derr2=cderr,id_in=0)
2279  call calctensord(td2,tduv2,tderr_aux2,cd,cduv,cderr,momvec,rmax)
2280 
2281  ! comparison --> take better result
2282  mode_cll = 3
2283  do r=0,rmax
2284  norm_coli=0d0
2285  norm_dd=0d0
2286  do n0=0,r
2287  do n1=0,r-n0
2288  do n2=0,r-n0-n1
2289  n3=r-n0-n1-n2
2290  norm_coli = max(norm_coli,abs(td(n0,n1,n2,n3)))
2291  norm_dd = max(norm_dd,abs(td2(n0,n1,n2,n3)))
2292  end do
2293  end do
2294  end do
2295  if (norm_coli.eq.0d0) then
2296  norm_coli = max(maxval(abs(mominv(1:6))),maxval(abs(masses2(0:3))))
2297  if(norm_coli.ne.0d0) then
2298  norm_coli=1d0/norm_coli**(2-real(r)/2)
2299  else
2300  norm_coli=1d0/muir2_cll**(2-real(r)/2)
2301  end if
2302  end if
2303  if (norm_dd.eq.0d0) then
2304  norm_dd = max(maxval(abs(mominv(1:6))),maxval(abs(masses2(0:3))))
2305  if(norm_dd.ne.0d0) then
2306  norm_dd=1d0/norm_dd**(2-real(r)/2)
2307  else
2308  norm_dd=1d0/muir2_cll**(2-real(r)/2)
2309  end if
2310  end if
2311  norm(r) = min(norm_coli,norm_dd)
2312  end do
2313 
2314  call checktensors_cll(td,td2,momvec,mominv,masses2,norm,4,rmax,tddiff)
2315 
2316  if (tderr_aux(rmax).lt.tderr_aux2(rmax)) then
2317  if (present(tderr)) tderr = max(tderr_aux,tddiff*norm)
2318  do r=0,rmax
2319  tdacc(r) = max(tderr_aux(r)/norm(r),tddiff(r))
2320  end do
2321  if (monitoring) pointscntdten_coli = pointscntdten_coli + 1
2322  else
2323  td = td2
2324  tduv = tduv2
2325  if (present(tderr)) tderr = max(tderr_aux2,tddiff*norm)
2326  do r=0,rmax
2327  tdacc(r) = max(tderr_aux2(r)/norm(r),tddiff(r))
2328  end do
2329  if (monitoring) pointscntdten_dd = pointscntdten_dd + 1
2330  end if
2331 
2332  else
2333  call d_main_cll(cd,cduv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
2334  masses2(0),masses2(1),masses2(2),masses2(3),rmax,derr2=cderr,id_in=0)
2335  call calctensord(td,tduv,tderr_aux,cd,cduv,cderr,momvec,rmax)
2336  if (present(tderr)) tderr = tderr_aux
2337  norm=0d0
2338  do r=0,rmax
2339  do n0=0,r
2340  do n1=0,r-n0
2341  do n2=0,r-n0-n1
2342  n3=r-n0-n1-n2
2343  norm(r) = max(norm(r),abs(td(n0,n1,n2,n3)))
2344  end do
2345  end do
2346  end do
2347  if (norm(r).eq.0d0) then
2348  norm(r) = max(maxval(abs(mominv(1:6))),maxval(abs(masses2(0:3))))
2349  if(norm(r).ne.0d0) then
2350  norm(r)=1d0/norm(r)**(2-real(r)/2)
2351  else
2352  norm(r)=1d0/muir2_cll**(2-real(r)/2)
2353  end if
2354  end if
2355  tdacc(r) = tderr_aux(r)/norm(r)
2356  end do
2357 
2358  end if
2359 
2360  call propagateaccflag_cll(tdacc,rmax)
2361  call propagateerrflag_cll
2362 
2363  if (monitoring) then
2364  pointscntdten_cll = pointscntdten_cll + 1
2365 
2366  if(maxval(tdacc).gt.reqacc_cll) accpointscntdten_cll = accpointscntdten_cll + 1
2367 
2368  if(maxval(tdacc).gt.critacc_cll) then
2369  critpointscntdten_cll = critpointscntdten_cll + 1
2370  if ( critpointscntdten_cll.le.noutcritpointsmax_cll(4) ) then
2371  call critpointsout_cll('TDten_cll',0,maxval(tdacc),critpointscntdten_cll)
2372  if( critpointscntdten_cll.eq.noutcritpointsmax_cll(4)) then
2373  write(ncpout_cll,*) ' Further output of Critical Points for TDten_cll suppressed'
2374  write(ncpout_cll,*)
2375  endif
2376 #ifdef CritPoints2
2377  call critpointsout2_cll('TDten_cll',0,maxval(tdacc),critpointscntdten_cll)
2378  if( critpointscntdten_cll.eq.noutcritpointsmax_cll(4)) then
2379  write(ncpout2_cll,*) ' Further output of Critical Points for TDten_cll suppressed'
2380  write(ncpout2_cll,*)
2381  endif
2382 #endif
2383  end if
2384  end if
2385  end if
2386 

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