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

Public Member Functions

subroutine eten_main_cll (TE, TEuv, MomVec, MomInv, masses2, rmax, TEerr)
 
subroutine eten_list_cll (TE, TEuv, MomVec, MomInv, masses2, rmax, TEerr)
 
subroutine eten_args_cll (TE, TEuv, p1vec, p2vec, p3vec, p4vec, p10, p21, p32, p43, p40, p20, p31, p42, p30, p41, m02, m12, m22, m32, m42, rmax, TEerr)
 
subroutine eten_args_list_cll (TE, TEuv, p1vec, p2vec, p3vec, p4vec, p10, p21, p32, p43, p40, p20, p31, p42, p30, p41, m02, m12, m22, m32, m42, rmax, TEerr)
 

Detailed Description

Definition at line 68 of file collier_tensors.F90.

Member Function/Subroutine Documentation

◆ eten_args_cll()

subroutine collier_tensors::eten_cll::eten_args_cll ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TE,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TEuv,
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, dimension(0:3), intent(in)  p4vec,
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  TEerr 
)

Definition at line 3341 of file collier_tensors.F90.

3341 
3342  integer, intent(in) :: rmax
3343  double complex, intent(in) :: p1vec(0:3),p2vec(0:3),p3vec(0:3),p4vec(0:3)
3344  double complex, intent(in) :: p10,p21,p32,p43,p40,p20,p31,p42,p30,p41
3345  double complex, intent(in) :: m02,m12,m22,m32,m42
3346  double complex, intent(out) :: TE(0:rmax,0:rmax,0:rmax,0:rmax)
3347  double complex, intent(out) :: TEuv(0:rmax,0:rmax,0:rmax,0:rmax)
3348  double precision, intent(out), optional :: TEerr(0:rmax)
3349  double complex :: TE2(0:rmax,0:rmax,0:rmax,0:rmax), TEuv2(0:rmax,0:rmax,0:rmax,0:rmax)
3350  double complex :: MomVec(0:3,4), MomInv(10), masses2(0:4)
3351  double complex :: CE(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
3352  double complex :: CEuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
3353  double precision :: CEerr(0:rmax), TEerr_aux(0:rmax), TEerr_aux2(0:rmax)
3354  double complex :: args(31)
3355  double precision :: TEdiff(0:rmax),norm(0:rmax),norm_coli,norm_dd,TEacc(0:rmax)
3356  integer :: r,n0,n1,n2,n3
3357  logical :: eflag
3358 
3359  if (5.gt.nmax_cll) then
3360  call seterrflag_cll(-10)
3361  call errout_cll('Eten_cll','Nmax_cll smaller 5',eflag,.true.)
3362  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
3363  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 5'
3364  call propagateerrflag_cll
3365  return
3366  end if
3367  if (rmax.gt.rmax_cll) then
3368  call seterrflag_cll(-10)
3369  call errout_cll('Eten_cll','argument rmax larger than rmax_cll',eflag,.true.)
3370  write(nerrout_cll,*) 'rmax =',rmax
3371  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
3372  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
3373  call propagateerrflag_cll
3374  return
3375  end if
3376 
3377  momvec(0:,1) = p1vec
3378  momvec(0:,2) = p2vec
3379  momvec(0:,3) = p3vec
3380  momvec(0:,4) = p4vec
3381  mominv(1) = p10
3382  mominv(2) = p21
3383  mominv(3) = p32
3384  mominv(4) = p43
3385  mominv(5) = p40
3386  mominv(6) = p20
3387  mominv(7) = p31
3388  mominv(8) = p42
3389  mominv(9) = p30
3390  mominv(10) = p41
3391  masses2(0) = m02
3392  masses2(1) = m12
3393  masses2(2) = m22
3394  masses2(3) = m32
3395  masses2(4) = m42
3396 
3397  ! set ID of master call
3398  args(1:4) = momvec(0:,1)
3399  args(5:8) = momvec(0:,2)
3400  args(9:12) = momvec(0:,3)
3401  args(13:16) = momvec(0:,4)
3402  args(17:26) = mominv
3403  args(27:31) = masses2
3404  call setmasterfname_cll('Eten_cll')
3405  call setmastern_cll(5)
3406  call setmasterr_cll(rmax)
3407  call setmasterargs_cll(31,args)
3408 
3409  call settencache_cll(tenred_cll-1)
3410 
3411 
3412  if (mode_cll.eq.3) then
3413  ! calculate tensor with coefficients from COLI
3414  mode_cll = 1
3415  call e_main_cll(ce,ceuv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5), &
3416  mominv(6),mominv(7),mominv(8),mominv(9),mominv(10),masses2(0), &
3417  masses2(1),masses2(2),masses2(3),masses2(4),rmax,eerr2=ceerr,id_in=0)
3418  call calctensore(te,teuv,teerr_aux,ce,ceuv,ceerr,momvec,rmax)
3419 
3420  ! calculate tensor with coefficients from DD
3421  mode_cll = 2
3422  call e_main_cll(ce,ceuv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5), &
3423  mominv(6),mominv(7),mominv(8),mominv(9),mominv(10),masses2(0), &
3424  masses2(1),masses2(2),masses2(3),masses2(4),rmax,eerr2=ceerr,id_in=0)
3425  call calctensore(te2,teuv2,teerr_aux2,ce,ceuv,ceerr,momvec,rmax)
3426 
3427  ! comparison --> take better result
3428  mode_cll = 3
3429  do r=0,rmax
3430  norm_coli=0d0
3431  norm_dd=0d0
3432  do n0=0,r
3433  do n1=0,r-n0
3434  do n2=0,r-n0-n1
3435  n3=r-n0-n1-n2
3436  norm_coli = max(norm_coli,abs(te(n0,n1,n2,n3)))
3437  norm_dd = max(norm_dd,abs(te2(n0,n1,n2,n3)))
3438  end do
3439  end do
3440  end do
3441  if (norm_coli.eq.0d0) then
3442  norm_coli = max(maxval(abs(mominv(1:10))),maxval(abs(masses2(0:4))))
3443  if(norm_coli.ne.0d0) then
3444  norm_coli=1d0/norm_coli**(3-real(r)/2)
3445  else
3446  norm_coli=1d0/muir2_cll**(3-real(r)/2)
3447  end if
3448  end if
3449  if (norm_dd.eq.0d0) then
3450  norm_dd = max(maxval(abs(mominv(1:10))),maxval(abs(masses2(0:4))))
3451  if(norm_dd.ne.0d0) then
3452  norm_dd=1d0/norm_dd**(3-real(r)/2)
3453  else
3454  norm_dd=1d0/muir2_cll**(3-real(r)/2)
3455  end if
3456  end if
3457  norm(r) = min(norm_coli,norm_dd)
3458  end do
3459 
3460  call checktensors_cll(te,te2,momvec,mominv,masses2,norm,5,rmax,tediff)
3461 
3462  if (teerr_aux(rmax).lt.teerr_aux2(rmax)) then
3463  if (present(teerr)) teerr = max(teerr_aux,tediff*norm)
3464  do r=0,rmax
3465  teacc(r) = max(teerr_aux(r)/norm(r),tediff(r))
3466  end do
3467  if (monitoring) pointscnteten_coli = pointscnteten_coli + 1
3468  else
3469  te = te2
3470  teuv = teuv2
3471  if (present(teerr)) teerr = max(teerr_aux2,tediff*norm)
3472  do r=0,rmax
3473  teacc(r) = max(teerr_aux2(r)/norm(r),tediff(r))
3474  end do
3475  if (monitoring) pointscnteten_dd = pointscnteten_dd + 1
3476  end if
3477 
3478  else
3479  call e_main_cll(ce,ceuv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5), &
3480  mominv(6),mominv(7),mominv(8),mominv(9),mominv(10),masses2(0), &
3481  masses2(1),masses2(2),masses2(3),masses2(4),rmax,eerr2=ceerr,id_in=0)
3482  call calctensore(te,teuv,teerr_aux,ce,ceuv,ceerr,momvec,rmax)
3483  if (present(teerr)) teerr = teerr_aux
3484  norm = 0d0
3485  do r=0,rmax
3486  do n0=0,r
3487  do n1=0,r-n0
3488  do n2=0,r-n0-n1
3489  n3=r-n0-n1-n2
3490  norm(r) = max(norm(r),abs(te(n0,n1,n2,n3)))
3491  end do
3492  end do
3493  end do
3494  if (norm(r).eq.0d0) then
3495  norm(r) = max(maxval(abs(mominv(1:10))),maxval(abs(masses2(0:4))))
3496  if(norm(r).ne.0d0) then
3497  norm(r)=1d0/norm(r)**(3-real(r)/2)
3498  else
3499  norm(r)=1d0/muir2_cll**(3-real(r)/2)
3500  end if
3501  end if
3502  teacc(r) = teerr_aux(r)/norm(r)
3503  end do
3504 
3505  end if
3506 
3507  call propagateaccflag_cll(teacc,rmax)
3508  call propagateerrflag_cll
3509 
3510  if (monitoring) then
3511  pointscnteten_cll = pointscnteten_cll + 1
3512 
3513  if(maxval(teacc).gt.reqacc_cll) accpointscnteten_cll = accpointscnteten_cll + 1
3514 
3515  if(maxval(teacc).gt.critacc_cll) then
3516  critpointscnteten_cll = critpointscnteten_cll + 1
3517  if ( critpointscnteten_cll.le.noutcritpointsmax_cll(5) ) then
3518  call critpointsout_cll('TEten_cll',0,maxval(teacc),critpointscnteten_cll)
3519  if( critpointscnteten_cll.eq.noutcritpointsmax_cll(5)) then
3520  write(ncpout_cll,*) ' Further output of Critical Points for TEten_cll suppressed'
3521  write(ncpout_cll,*)
3522  endif
3523 #ifdef CritPoints2
3524  call critpointsout2_cll('TEten_cll',0,maxval(teacc),critpointscnteten_cll)
3525  if( critpointscnteten_cll.eq.noutcritpointsmax_cll(5)) then
3526  write(ncpout2_cll,*) ' Further output of Critical Points for TEten_cll suppressed'
3527  write(ncpout2_cll,*)
3528  endif
3529 #endif
3530  end if
3531  end if
3532  end if
3533 

◆ eten_args_list_cll()

subroutine collier_tensors::eten_cll::eten_args_list_cll ( double complex, dimension(rts(rmax)), intent(out)  TE,
double complex, dimension(rts(rmax)), intent(out)  TEuv,
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, dimension(0:3), intent(in)  p4vec,
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  TEerr 
)

Definition at line 3548 of file collier_tensors.F90.

3548 
3549  integer, intent(in) :: rmax
3550  double complex, intent(in) :: p1vec(0:3),p2vec(0:3),p3vec(0:3),p4vec(0:3)
3551  double complex, intent(in) :: p10,p21,p32,p43,p40,p20,p31,p42,p30,p41
3552  double complex, intent(in) :: m02,m12,m22,m32,m42
3553  double complex, intent(out) :: TE(RtS(rmax)), TEuv(RtS(rmax))
3554  double precision, intent(out), optional :: TEerr(0:rmax)
3555  double complex :: TE2(RtS(rmax)), TEuv2(RtS(rmax))
3556  double complex :: MomVec(0:3,4), MomInv(10), masses2(0:4)
3557  double complex :: CE(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
3558  double complex :: CEuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
3559  double precision :: CEerr(0:rmax), TEerr_aux(0:rmax), TEerr_aux2(0:rmax)
3560  double complex :: args(31)
3561  double precision :: TEdiff(0:rmax),norm(0:rmax),norm_coli,norm_dd,TEacc(0:rmax)
3562  integer :: r,i
3563  logical :: eflag
3564 
3565  if (5.gt.nmax_cll) then
3566  call seterrflag_cll(-10)
3567  call errout_cll('Eten_cll','Nmax_cll smaller 5',eflag,.true.)
3568  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
3569  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 5'
3570  call propagateerrflag_cll
3571  return
3572  end if
3573  if (rmax.gt.rmax_cll) then
3574  call seterrflag_cll(-10)
3575  call errout_cll('Eten_cll','argument rmax larger than rmax_cll',eflag,.true.)
3576  write(nerrout_cll,*) 'rmax =',rmax
3577  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
3578  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
3579  call propagateerrflag_cll
3580  return
3581  end if
3582 
3583  call eten_args_list_checked_cll(te,teuv,p1vec,p2vec,p3vec,p4vec, &
3584  p10,p21,p32,p43,p40,p20,p31,p42,p30,p41, &
3585  m02,m12,m22,m32,m42,rmax,teerr)
3586 

◆ eten_list_cll()

subroutine collier_tensors::eten_cll::eten_list_cll ( double complex, dimension(:), intent(out)  TE,
double complex, dimension(:), intent(out)  TEuv,
double complex, dimension(0:3,4), intent(in)  MomVec,
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  TEerr 
)

Definition at line 3155 of file collier_tensors.F90.

3155 
3156  integer, intent(in) :: rmax
3157  double complex, intent(in) :: MomVec(0:3,4), MomInv(10), masses2(0:4)
3158  double complex, intent(out) :: TE(:), TEuv(:)
3159  double precision, intent(out), optional :: TEerr(0:rmax)
3160  integer :: r,i
3161  logical :: eflag
3162 
3163  if (5.gt.nmax_cll) then
3164  call seterrflag_cll(-10)
3165  call errout_cll('Eten_cll','Nmax_cll smaller 5',eflag,.true.)
3166  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
3167  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 5'
3168  call propagateerrflag_cll
3169  return
3170  end if
3171  if (rmax.gt.rmax_cll) then
3172  call seterrflag_cll(-10)
3173  call errout_cll('Eten_cll','argument rmax larger than rmax_cll',eflag,.true.)
3174  write(nerrout_cll,*) 'rmax =',rmax
3175  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
3176  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
3177  call propagateerrflag_cll
3178  return
3179  end if
3180 
3181  call eten_list_checked_cll(te,teuv,momvec,mominv,masses2,rmax,teerr)
3182 

◆ eten_main_cll()

subroutine collier_tensors::eten_cll::eten_main_cll ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TE,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TEuv,
double complex, dimension(0:3,4), intent(in)  MomVec,
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  TEerr 
)

Definition at line 2974 of file collier_tensors.F90.

2974 
2975  integer, intent(in) :: rmax
2976  double complex, intent(in) :: MomVec(0:3,4), MomInv(10), masses2(0:4)
2977  double complex, intent(out) :: TE(0:rmax,0:rmax,0:rmax,0:rmax)
2978  double complex, intent(out) :: TEuv(0:rmax,0:rmax,0:rmax,0:rmax)
2979  double precision, intent(out), optional :: TEerr(0:rmax)
2980  double complex :: TE2(0:rmax,0:rmax,0:rmax,0:rmax), TEuv2(0:rmax,0:rmax,0:rmax,0:rmax)
2981  double complex :: CE(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
2982  double complex :: CEuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax)
2983  double precision :: CEerr(0:rmax), TEerr_aux(0:rmax), TEerr_aux2(0:rmax)
2984  double complex :: args(31)
2985  double precision :: TEdiff(0:rmax),norm(0:rmax),norm_coli,norm_dd,TEacc(0:rmax)
2986  integer :: r,n0,n1,n2,n3
2987  logical :: eflag
2988 
2989  if (5.gt.nmax_cll) then
2990  call seterrflag_cll(-10)
2991  call errout_cll('Eten_cll','Nmax_cll smaller 5',eflag,.true.)
2992  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
2993  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 5'
2994  call propagateerrflag_cll
2995  return
2996  end if
2997  if (rmax.gt.rmax_cll) then
2998  call seterrflag_cll(-10)
2999  call errout_cll('Eten_cll','argument rmax larger than rmax_cll',eflag,.true.)
3000  write(nerrout_cll,*) 'rmax =',rmax
3001  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
3002  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
3003  call propagateerrflag_cll
3004  return
3005  end if
3006 
3007  ! set ID of master call
3008  args(1:4) = momvec(0:,1)
3009  args(5:8) = momvec(0:,2)
3010  args(9:12) = momvec(0:,3)
3011  args(13:16) = momvec(0:,4)
3012  args(17:26) = mominv
3013  args(27:31) = masses2
3014  call setmasterfname_cll('Eten_cll')
3015  call setmastern_cll(5)
3016  call setmasterr_cll(rmax)
3017  call setmasterargs_cll(31,args)
3018 
3019  call settencache_cll(tenred_cll-1)
3020 
3021  if (mode_cll.eq.3) then
3022  ! calculate tensor with coefficients from COLI
3023  mode_cll = 1
3024  call e_main_cll(ce,ceuv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5), &
3025  mominv(6),mominv(7),mominv(8),mominv(9),mominv(10),masses2(0), &
3026  masses2(1),masses2(2),masses2(3),masses2(4),rmax,eerr2=ceerr,id_in=0)
3027  call calctensore(te,teuv,teerr_aux,ce,ceuv,ceerr,momvec,rmax)
3028 
3029  ! calculate tensor with coefficients from DD
3030  mode_cll = 2
3031  call e_main_cll(ce,ceuv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5), &
3032  mominv(6),mominv(7),mominv(8),mominv(9),mominv(10),masses2(0), &
3033  masses2(1),masses2(2),masses2(3),masses2(4),rmax,eerr2=ceerr,id_in=0)
3034  call calctensore(te2,teuv2,teerr_aux2,ce,ceuv,ceerr,momvec,rmax)
3035 
3036  ! comparison --> take better result
3037  mode_cll = 3
3038  do r=0,rmax
3039  norm_coli=0d0
3040  norm_dd=0d0
3041  do n0=0,r
3042  do n1=0,r-n0
3043  do n2=0,r-n0-n1
3044  n3=r-n0-n1-n2
3045  norm_coli = max(norm_coli,abs(te(n0,n1,n2,n3)))
3046  norm_dd = max(norm_dd,abs(te2(n0,n1,n2,n3)))
3047  end do
3048  end do
3049  end do
3050  if (norm_coli.eq.0d0) then
3051  norm_coli = max(maxval(abs(mominv(1:10))),maxval(abs(masses2(0:4))))
3052  if(norm_coli.ne.0d0) then
3053  norm_coli=1d0/norm_coli**(3-real(r)/2)
3054  else
3055  norm_coli=1d0/muir2_cll**(3-real(r)/2)
3056  end if
3057  end if
3058  if (norm_dd.eq.0d0) then
3059  norm_dd = max(maxval(abs(mominv(1:10))),maxval(abs(masses2(0:4))))
3060  if(norm_dd.ne.0d0) then
3061  norm_dd=1d0/norm_dd**(3-real(r)/2)
3062  else
3063  norm_dd=1d0/muir2_cll**(3-real(r)/2)
3064  end if
3065  end if
3066  norm(r) = min(norm_coli,norm_dd)
3067  end do
3068 
3069  call checktensors_cll(te,te2,momvec,mominv,masses2,norm,5,rmax,tediff)
3070 
3071  if (teerr_aux(rmax).lt.teerr_aux2(rmax)) then
3072  if (present(teerr)) teerr = max(teerr_aux,tediff*norm)
3073  do r=0,rmax
3074  teacc(r) = max(teerr_aux(r)/norm(r),tediff(r))
3075  end do
3076  if (monitoring) pointscnteten_coli = pointscnteten_coli + 1
3077  else
3078  te = te2
3079  teuv = teuv2
3080  if (present(teerr)) teerr = max(teerr_aux2,tediff*norm)
3081  do r=0,rmax
3082  teacc(r) = max(teerr_aux2(r)/norm(r),tediff(r))
3083  end do
3084  if (monitoring) pointscnteten_dd = pointscnteten_dd + 1
3085  end if
3086 
3087  else
3088  call e_main_cll(ce,ceuv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5), &
3089  mominv(6),mominv(7),mominv(8),mominv(9),mominv(10),masses2(0), &
3090  masses2(1),masses2(2),masses2(3),masses2(4),rmax,eerr2=ceerr,id_in=0)
3091  call calctensore(te,teuv,teerr_aux,ce,ceuv,ceerr,momvec,rmax)
3092  if (present(teerr)) teerr = teerr_aux
3093  norm = 0d0
3094  do r=0,rmax
3095  do n0=0,r
3096  do n1=0,r-n0
3097  do n2=0,r-n0-n1
3098  n3=r-n0-n1-n2
3099  norm(r) = max(norm(r),abs(te(n0,n1,n2,n3)))
3100  end do
3101  end do
3102  end do
3103  if (norm(r).eq.0d0) then
3104  norm(r) = max(maxval(abs(mominv(1:10))),maxval(abs(masses2(0:4))))
3105  if(norm(r).ne.0d0) then
3106  norm(r)=1d0/norm(r)**(3-real(r)/2)
3107  else
3108  norm(r)=1d0/muir2_cll**(3-real(r)/2)
3109  end if
3110  end if
3111  teacc(r) = teerr_aux(r)/norm(r)
3112  end do
3113 
3114  end if
3115 
3116  call propagateaccflag_cll(teacc,rmax)
3117  call propagateerrflag_cll
3118 
3119  if (monitoring) then
3120  pointscnteten_cll = pointscnteten_cll + 1
3121 
3122  if(maxval(teacc).gt.reqacc_cll) accpointscnteten_cll = accpointscnteten_cll + 1
3123 
3124  if(maxval(teacc).gt.critacc_cll) then
3125  critpointscnteten_cll = critpointscnteten_cll + 1
3126  if ( critpointscnteten_cll.le.noutcritpointsmax_cll(5) ) then
3127  call critpointsout_cll('TEten_cll',0,maxval(teacc),critpointscnteten_cll)
3128  if( critpointscnteten_cll.eq.noutcritpointsmax_cll(5)) then
3129  write(ncpout_cll,*) ' Further output of Critical Points for TEten_cll suppressed'
3130  write(ncpout_cll,*)
3131  endif
3132 #ifdef CritPoints2
3133  call critpointsout2_cll('TEten_cll',0,maxval(teacc),critpointscnteten_cll)
3134  if( critpointscnteten_cll.eq.noutcritpointsmax_cll(5)) then
3135  write(ncpout2_cll,*) ' Further output of Critical Points for TEten_cll suppressed'
3136  write(ncpout2_cll,*)
3137  endif
3138 #endif
3139  end if
3140  end if
3141  end if
3142 

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