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

Public Member Functions

subroutine gten_main_cll (TG, TGuv, MomVec, MomInv, masses2, rmax, TGerr)
 
subroutine gten_list_cll (TG, TGuv, MomVec, MomInv, masses2, rmax, TGerr)
 
subroutine gten_args_cll (TG, TGuv, p1vec, p2vec, p3vec, p4vec, p5vec, p6vec, p10, p21, p32, p43, p54, p65, p60, p20, p31, p42, p53, p64, p50, p61, p30, p41, p52, p63, p40, p51, p62, m02, m12, m22, m32, m42, m52, m62, rmax, TGerr)
 
subroutine gten_args_list_cll (TG, TGuv, p1vec, p2vec, p3vec, p4vec, p5vec, p6vec, p10, p21, p32, p43, p54, p65, p60, p20, p31, p42, p53, p64, p50, p61, p30, p41, p52, p63, p40, p51, p62, m02, m12, m22, m32, m42, m52, m62, rmax, TGerr)
 

Detailed Description

Definition at line 80 of file collier_tensors.F90.

Member Function/Subroutine Documentation

◆ gten_args_cll()

subroutine collier_tensors::gten_cll::gten_args_cll ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TG,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TGuv,
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, dimension(0:3), intent(in)  p5vec,
double complex, dimension(0:3), intent(in)  p6vec,
double complex, intent(in)  p10,
double complex, intent(in)  p21,
double complex, intent(in)  p32,
double complex, intent(in)  p43,
double complex, intent(in)  p54,
double complex, intent(in)  p65,
double complex, intent(in)  p60,
double complex, intent(in)  p20,
double complex, intent(in)  p31,
double complex, intent(in)  p42,
double complex, intent(in)  p53,
double complex, intent(in)  p64,
double complex, intent(in)  p50,
double complex, intent(in)  p61,
double complex, intent(in)  p30,
double complex, intent(in)  p41,
double complex, intent(in)  p52,
double complex, intent(in)  p63,
double complex, intent(in)  p40,
double complex, intent(in)  p51,
double complex, intent(in)  p62,
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,
double complex, intent(in)  m52,
double complex, intent(in)  m62,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  TGerr 
)

Definition at line 5371 of file collier_tensors.F90.

5371  integer, intent(in) :: rmax
5372  double complex, intent(in) :: p1vec(0:3),p2vec(0:3),p3vec(0:3),p4vec(0:3)
5373  double complex, intent(in) :: p5vec(0:3),p6vec(0:3)
5374  double complex, intent(in) :: p10,p21,p32,p43,p54,p65,p60,p20,p31,p42,p53
5375  double complex, intent(in) :: p64,p50,p61,p30,p41,p52,p63,p40,p51,p62
5376  double complex, intent(in) :: m02,m12,m22,m32,m42,m52,m62
5377  double complex, intent(out) :: TG(0:rmax,0:rmax,0:rmax,0:rmax)
5378  double complex, intent(out) :: TGuv(0:rmax,0:rmax,0:rmax,0:rmax)
5379  double precision, intent(out), optional :: TGerr(0:rmax)
5380  double complex :: TG2(0:rmax,0:rmax,0:rmax,0:rmax), TGuv2(0:rmax,0:rmax,0:rmax,0:rmax)
5381  double precision :: TGerr_aux(0:rmax),TGerr_aux2(0:rmax)
5382  double complex :: MomVec(0:3,6), MomInv(21), masses2(0:6)
5383  double complex :: CG(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
5384  double complex :: CGuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
5385  double precision :: CGerr(0:rmax), TGacc(0:rmax)
5386  double precision :: norm(0:rmax),norm_coli,norm_dd, TGdiff(0:rmax)
5387  double complex :: args(52)
5388  integer :: r,n0,n1,n2,n3
5389  logical :: eflag
5390 
5391  if (7.gt.nmax_cll) then
5392  call seterrflag_cll(-10)
5393  call errout_cll('Gten_cll','Nmax_cll smaller 7',eflag,.true.)
5394  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
5395  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 7'
5396  call propagateerrflag_cll
5397  return
5398  end if
5399  if (rmax.gt.rmax_cll) then
5400  call seterrflag_cll(-10)
5401  call errout_cll('Gten_cll','argument rmax larger than rmax_cll',eflag,.true.)
5402  write(nerrout_cll,*) 'rmax =',rmax
5403  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
5404  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
5405  call propagateerrflag_cll
5406  return
5407  end if
5408 
5409  momvec(0:,1) = p1vec
5410  momvec(0:,2) = p2vec
5411  momvec(0:,3) = p3vec
5412  momvec(0:,4) = p4vec
5413  momvec(0:,5) = p5vec
5414  momvec(0:,6) = p6vec
5415  mominv(1) = p10
5416  mominv(2) = p21
5417  mominv(3) = p32
5418  mominv(4) = p43
5419  mominv(5) = p54
5420  mominv(6) = p65
5421  mominv(7) = p60
5422  mominv(8) = p20
5423  mominv(9) = p31
5424  mominv(10) = p42
5425  mominv(11) = p53
5426  mominv(12) = p64
5427  mominv(13) = p50
5428  mominv(14) = p61
5429  mominv(15) = p30
5430  mominv(16) = p41
5431  mominv(17) = p52
5432  mominv(18) = p63
5433  mominv(19) = p40
5434  mominv(20) = p51
5435  mominv(21) = p62
5436  masses2(0) = m02
5437  masses2(1) = m12
5438  masses2(2) = m22
5439  masses2(3) = m32
5440  masses2(4) = m42
5441  masses2(5) = m52
5442  masses2(6) = m62
5443 
5444  ! set ID of master call
5445  args(1:4) = momvec(0:,1)
5446  args(5:8) = momvec(0:,2)
5447  args(9:12) = momvec(0:,3)
5448  args(13:16) = momvec(0:,4)
5449  args(17:20) = momvec(0:,5)
5450  args(21:24) = momvec(0:,6)
5451  args(25:45) = mominv
5452  args(46:52) = masses2
5453  call setmasterfname_cll('Gten_cll')
5454  call setmastern_cll(7)
5455  call setmasterr_cll(rmax)
5456  call setmasterargs_cll(52,args)
5457 
5458  call settencache_cll(tenred_cll-1)
5459 
5460 
5461  if (tenred_cll.le.7) then
5462 
5463  if (mode_cll.gt.1) call tn_dd_dummy(7,rmax)
5464 
5465  if (mode_cll.eq.3) then
5466  ! calculate tensor with coefficients from COLI
5467  mode_cll = 1
5468  call calctensortnr(tg,tguv,tgerr_aux,momvec,mominv,masses2,7,rmax,0)
5469 
5470  ! calculate tensor with coefficients from DD
5471  mode_cll = 2
5472  call calctensortnr(tg2,tguv2,tgerr_aux2,momvec,mominv,masses2,7,rmax,0)
5473 
5474  ! comparison --> take better result
5475  mode_cll = 3
5476  do r=0,rmax
5477  norm_coli=0d0
5478  norm_dd=0d0
5479  do n0=0,r
5480  do n1=0,r-n0
5481  do n2=0,r-n0-n1
5482  n3=r-n0-n1-n2
5483  norm_coli = max(norm_coli,abs(tg(n0,n1,n2,n3)))
5484  norm_dd = max(norm_dd,abs(tg2(n0,n1,n2,n3)))
5485  end do
5486  end do
5487  end do
5488  if (norm_coli.eq.0d0) then
5489  norm_coli = max(maxval(abs(mominv(1:21))),maxval(abs(masses2(0:6))))
5490  if(norm_coli.ne.0d0) then
5491  norm_coli=1d0/norm_coli**(5-real(r)/2)
5492  else
5493  norm_coli=1d0/muir2_cll**(5-real(r)/2)
5494  end if
5495  end if
5496  if (norm_dd.eq.0d0) then
5497  norm_dd = max(maxval(abs(mominv(1:21))),maxval(abs(masses2(0:6))))
5498  if(norm_dd.ne.0d0) then
5499  norm_dd=1d0/norm_dd**(5-real(r)/2)
5500  else
5501  norm_dd=1d0/muir2_cll**(5-real(r)/2)
5502  end if
5503  end if
5504  norm(r) = min(norm_coli,norm_dd)
5505  end do
5506 
5507  call checktensors_cll(tg,tg2,momvec,mominv,masses2,norm,7,rmax,tgdiff)
5508 
5509  if (tgerr_aux(rmax).lt.tgerr_aux2(rmax)) then
5510  if (present(tgerr)) tgerr = max(tgerr_aux,tgdiff*norm)
5511  do r=0,rmax
5512  tgacc(r) = max(tgerr_aux(r)/norm(r),tgdiff(r))
5513  end do
5514  if (monitoring) pointscntgten_coli = pointscntgten_coli + 1
5515  else
5516  tg = tg2
5517  tguv = tguv2
5518  if (present(tgerr)) tgerr = max(tgerr_aux2,tgdiff*norm)
5519  do r=0,rmax
5520  tgacc(r) = max(tgerr_aux2(r)/norm(r),tgdiff(r))
5521  end do
5522  if (monitoring) pointscntgten_dd = pointscntgten_dd + 1
5523  end if
5524 
5525  else
5526  call calctensortnr(tg,tguv,tgerr_aux,momvec,mominv,masses2,7,rmax,0)
5527  if (present(tgerr)) tgerr = tgerr_aux
5528  norm = 0d0
5529  do r=0,rmax
5530  do n0=0,r
5531  do n1=0,r-n0
5532  do n2=0,r-n0-n1
5533  n3=r-n0-n1-n2
5534  norm(r) = max(norm(r),abs(tg(n0,n1,n2,n3)))
5535  end do
5536  end do
5537  end do
5538  if (norm(r).eq.0d0) then
5539  norm(r) = max(maxval(abs(mominv(1:21))),maxval(abs(masses2(0:6))))
5540  if(norm(r).ne.0d0) then
5541  norm(r)=1d0/norm(r)**(5-real(r)/2)
5542  else
5543  norm(r)=1d0/muir2_cll**(5-real(r)/2)
5544  end if
5545  end if
5546  tgacc(r) = tgerr_aux(r)/norm(r)
5547  end do
5548 
5549  end if
5550 
5551  else
5552  call g_main_cll(cg,cguv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
5553  mominv(7),mominv(8),mominv(9),mominv(10),mominv(11),mominv(12), &
5554  mominv(13),mominv(14),mominv(15),mominv(16),mominv(17),mominv(18), &
5555  mominv(19),mominv(20),mominv(21),masses2(0),masses2(1),masses2(2), &
5556  masses2(3),masses2(4),masses2(5),masses2(6),rmax,gerr2=cgerr,id_in=0)
5557  call calctensorg(tg,tguv,tgerr_aux,cg,cguv,cgerr,momvec,rmax)
5558  if (present(tgerr)) tgerr = tgerr_aux
5559  norm = 0d0
5560  do r=0,rmax
5561  do n0=0,r
5562  do n1=0,r-n0
5563  do n2=0,r-n0-n1
5564  n3=r-n0-n1-n2
5565  norm(r) = max(norm(r),abs(tg(n0,n1,n2,n3)))
5566  end do
5567  end do
5568  end do
5569  if (norm(r).eq.0d0) then
5570  norm(r) = max(maxval(abs(mominv(1:21))),maxval(abs(masses2(0:6))))
5571  if(norm(r).ne.0d0) then
5572  norm(r)=1d0/norm(r)**(5-real(r)/2)
5573  else
5574  norm(r)=1d0/muir2_cll**(5-real(r)/2)
5575  end if
5576  end if
5577  tgacc(r) = tgerr_aux(r)/norm(r)
5578  end do
5579  end if
5580 
5581  if (monitoring) then
5582  pointscntgten_cll = pointscntgten_cll + 1
5583 
5584  if(maxval(tgacc).gt.reqacc_cll) accpointscntgten_cll = accpointscntgten_cll + 1
5585 
5586  if(maxval(tgacc).gt.critacc_cll) then
5587  critpointscntgten_cll = critpointscntgten_cll + 1
5588  if ( critpointscntgten_cll.le.noutcritpointsmax_cll(7) ) then
5589  call critpointsout_cll('TGten_cll',0,maxval(tgacc),critpointscntgten_cll)
5590  if( critpointscntgten_cll.eq.noutcritpointsmax_cll(7)) then
5591  write(ncpout_cll,*) ' Further output of Critical Points for TGten_cll suppressed'
5592  write(ncpout_cll,*)
5593  endif
5594  end if
5595  end if
5596  end if
5597 

◆ gten_args_list_cll()

subroutine collier_tensors::gten_cll::gten_args_list_cll ( double complex, dimension(rts(rmax)), intent(out)  TG,
double complex, dimension(rts(rmax)), intent(out)  TGuv,
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, dimension(0:3), intent(in)  p5vec,
double complex, dimension(0:3), intent(in)  p6vec,
double complex, intent(in)  p10,
double complex, intent(in)  p21,
double complex, intent(in)  p32,
double complex, intent(in)  p43,
double complex, intent(in)  p54,
double complex, intent(in)  p65,
double complex, intent(in)  p60,
double complex, intent(in)  p20,
double complex, intent(in)  p31,
double complex, intent(in)  p42,
double complex, intent(in)  p53,
double complex, intent(in)  p64,
double complex, intent(in)  p50,
double complex, intent(in)  p61,
double complex, intent(in)  p30,
double complex, intent(in)  p41,
double complex, intent(in)  p52,
double complex, intent(in)  p63,
double complex, intent(in)  p40,
double complex, intent(in)  p51,
double complex, intent(in)  p62,
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,
double complex, intent(in)  m52,
double complex, intent(in)  m62,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  TGerr 
)

Definition at line 5616 of file collier_tensors.F90.

5616  integer, intent(in) :: rmax
5617  double complex, intent(in) :: p1vec(0:3),p2vec(0:3),p3vec(0:3),p4vec(0:3)
5618  double complex, intent(in) :: p5vec(0:3),p6vec(0:3)
5619  double complex, intent(in) :: p10,p21,p32,p43,p54,p65,p60,p20,p31,p42,p53
5620  double complex, intent(in) :: p64,p50,p61,p30,p41,p52,p63,p40,p51,p62
5621  double complex, intent(in) :: m02,m12,m22,m32,m42,m52,m62
5622  double complex, intent(out) :: TG(RtS(rmax)), TGuv(RtS(rmax))
5623  double precision, intent(out), optional :: TGerr(0:rmax)
5624  double complex :: TG2(RtS(rmax)), TGuv2(RtS(rmax))
5625  double precision :: TGerr_aux(0:rmax), TGerr_aux2(0:rmax)
5626  double complex :: MomVec(0:3,6), MomInv(21), masses2(0:6)
5627  double complex :: CG(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
5628  double complex :: CGuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
5629  double precision :: CGerr(0:rmax), TGacc(0:rmax)
5630  double precision :: norm(0:rmax), TGdiff(0:rmax), norm_coli, norm_dd
5631  double complex :: args(52)
5632  integer :: r,i
5633  logical :: eflag
5634 
5635  if (7.gt.nmax_cll) then
5636  call seterrflag_cll(-10)
5637  call errout_cll('Gten_cll','Nmax_cll smaller 7',eflag,.true.)
5638  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
5639  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 7'
5640  call propagateerrflag_cll
5641  return
5642  end if
5643  if (rmax.gt.rmax_cll) then
5644  call seterrflag_cll(-10)
5645  call errout_cll('Gten_cll','argument rmax larger than rmax_cll',eflag,.true.)
5646  write(nerrout_cll,*) 'rmax =',rmax
5647  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
5648  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
5649  call propagateerrflag_cll
5650  return
5651  end if
5652 
5653  call gten_args_list_checked_cll(tg,tguv,p1vec,p2vec,p3vec,p4vec,p5vec,p6vec, &
5654  p10,p21,p32,p43,p54,p65,p60,p20,p31,p42,p53, &
5655  p64,p50,p61,p30,p41,p52,p63,p40,p51,p62, &
5656  m02,m12,m22,m32,m42,m52,m62,rmax,tgerr)
5657 

◆ gten_list_cll()

subroutine collier_tensors::gten_cll::gten_list_cll ( double complex, dimension(:), intent(out)  TG,
double complex, dimension(:), intent(out)  TGuv,
double complex, dimension(0:3,6), intent(in)  MomVec,
double complex, dimension(21), intent(in)  MomInv,
double complex, dimension(0:6), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  TGerr 
)

Definition at line 5168 of file collier_tensors.F90.

5168 
5169  integer, intent(in) :: rmax
5170  double complex, intent(in) :: MomVec(0:3,6), MomInv(21), masses2(0:6)
5171  double complex, intent(out) :: TG(:),TGuv(:)
5172  double precision, intent(out), optional :: TGerr(0:rmax)
5173  logical :: eflag
5174 
5175  if (7.gt.nmax_cll) then
5176  call seterrflag_cll(-10)
5177  call errout_cll('Gten_cll','Nmax_cll smaller 7',eflag,.true.)
5178  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
5179  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 7'
5180  call propagateerrflag_cll
5181  return
5182  end if
5183  if (rmax.gt.rmax_cll) then
5184  call seterrflag_cll(-10)
5185  call errout_cll('Gten_cll','argument rmax larger than rmax_cll',eflag,.true.)
5186  write(nerrout_cll,*) 'rmax =',rmax
5187  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
5188  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
5189  call propagateerrflag_cll
5190  return
5191  end if
5192 
5193  call gten_list_checked_cll(tg,tguv,momvec,mominv,masses2,rmax,tgerr)
5194 

◆ gten_main_cll()

subroutine collier_tensors::gten_cll::gten_main_cll ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TG,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TGuv,
double complex, dimension(0:3,6), intent(in)  MomVec,
double complex, dimension(21), intent(in)  MomInv,
double complex, dimension(0:6), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  TGerr 
)

Definition at line 4945 of file collier_tensors.F90.

4945 
4946  integer, intent(in) :: rmax
4947  double complex, intent(in) :: MomVec(0:3,6), MomInv(21), masses2(0:6)
4948  double complex, intent(out) :: TG(0:rmax,0:rmax,0:rmax,0:rmax)
4949  double complex, intent(out) :: TGuv(0:rmax,0:rmax,0:rmax,0:rmax)
4950  double precision, intent(out), optional :: TGerr(0:rmax)
4951  double precision :: TGerr_aux(0:rmax), TGerr_aux2(0:rmax)
4952  double complex :: TG2(0:rmax,0:rmax,0:rmax,0:rmax), TGuv2(0:rmax,0:rmax,0:rmax,0:rmax)
4953  double complex :: CG(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
4954  double complex :: CGuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
4955  double precision :: CGerr(0:rmax), TGacc(0:rmax)
4956  double precision :: norm(0:rmax),norm_coli,norm_dd, TGdiff(0:rmax)
4957  double complex :: args(52)
4958  integer :: r,n0,n1,n2,n3
4959  logical :: eflag
4960 
4961  if (7.gt.nmax_cll) then
4962  call seterrflag_cll(-10)
4963  call errout_cll('Gten_cll','Nmax_cll smaller 7',eflag,.true.)
4964  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
4965  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 7'
4966  call propagateerrflag_cll
4967  return
4968  end if
4969  if (rmax.gt.rmax_cll) then
4970  call seterrflag_cll(-10)
4971  call errout_cll('Gten_cll','argument rmax larger than rmax_cll',eflag,.true.)
4972  write(nerrout_cll,*) 'rmax =',rmax
4973  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
4974  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
4975  call propagateerrflag_cll
4976  return
4977  end if
4978 
4979  ! set ID of master call
4980  args(1:4) = momvec(0:,1)
4981  args(5:8) = momvec(0:,2)
4982  args(9:12) = momvec(0:,3)
4983  args(13:16) = momvec(0:,4)
4984  args(17:20) = momvec(0:,5)
4985  args(21:24) = momvec(0:,6)
4986  args(25:45) = mominv
4987  args(46:52) = masses2
4988  call setmasterfname_cll('Gten_cll')
4989  call setmastern_cll(7)
4990  call setmasterr_cll(rmax)
4991  call setmasterargs_cll(52,args)
4992 
4993  call settencache_cll(tenred_cll-1)
4994 
4995  if (tenred_cll.le.7) then
4996 
4997  if (mode_cll.gt.1) call tn_dd_dummy(7,rmax)
4998 
4999  if (mode_cll.eq.3) then
5000  ! calculate tensor with coefficients from COLI
5001  mode_cll = 1
5002  call calctensortnr(tg,tguv,tgerr_aux,momvec,mominv,masses2,7,rmax,0)
5003 
5004  ! calculate tensor with coefficients from DD
5005  mode_cll = 2
5006  call calctensortnr(tg2,tguv2,tgerr_aux2,momvec,mominv,masses2,7,rmax,0)
5007 
5008  ! comparison --> take better result
5009  mode_cll = 3
5010  do r=0,rmax
5011  norm_coli=0d0
5012  norm_dd=0d0
5013  do n0=0,r
5014  do n1=0,r-n0
5015  do n2=0,r-n0-n1
5016  n3=r-n0-n1-n2
5017  norm_coli = max(norm_coli,abs(tg(n0,n1,n2,n3)))
5018  norm_dd = max(norm_dd,abs(tg2(n0,n1,n2,n3)))
5019  end do
5020  end do
5021  end do
5022  if (norm_coli.eq.0d0) then
5023  norm_coli = max(maxval(abs(mominv(1:21))),maxval(abs(masses2(0:6))))
5024  if(norm_coli.ne.0d0) then
5025  norm_coli=1d0/norm_coli**(5-real(r)/2)
5026  else
5027  norm_coli=1d0/muir2_cll**(5-real(r)/2)
5028  end if
5029  end if
5030  if (norm_dd.eq.0d0) then
5031  norm_dd = max(maxval(abs(mominv(1:21))),maxval(abs(masses2(0:6))))
5032  if(norm_dd.ne.0d0) then
5033  norm_dd=1d0/norm_dd**(5-real(r)/2)
5034  else
5035  norm_dd=1d0/muir2_cll**(5-real(r)/2)
5036  end if
5037  end if
5038  norm(r) = min(norm_coli,norm_dd)
5039  end do
5040 
5041  call checktensors_cll(tg,tg2,momvec,mominv,masses2,norm,7,rmax,tgdiff)
5042 
5043  if (tgerr_aux(rmax).lt.tgerr_aux2(rmax)) then
5044  if (present(tgerr)) tgerr = max(tgerr_aux,tgdiff*norm)
5045  do r=0,rmax
5046  tgacc(r) = max(tgerr_aux(r)/norm(r),tgdiff(r))
5047  end do
5048  if (monitoring) pointscntgten_coli = pointscntgten_coli + 1
5049  else
5050  tg = tg2
5051  tguv = tguv2
5052  if (present(tgerr)) tgerr = max(tgerr_aux2,tgdiff*norm)
5053  do r=0,rmax
5054  tgacc(r) = max(tgerr_aux2(r)/norm(r),tgdiff(r))
5055  end do
5056  if (monitoring) pointscntgten_dd = pointscntgten_dd + 1
5057  end if
5058 
5059  else
5060  call calctensortnr(tg,tguv,tgerr_aux,momvec,mominv,masses2,7,rmax,0)
5061  if (present(tgerr)) tgerr = tgerr_aux
5062  norm = 0d0
5063  do r=0,rmax
5064  do n0=0,r
5065  do n1=0,r-n0
5066  do n2=0,r-n0-n1
5067  n3=r-n0-n1-n2
5068  norm(r) = max(norm(r),abs(tg(n0,n1,n2,n3)))
5069  end do
5070  end do
5071  end do
5072  if (norm(r).eq.0d0) then
5073  norm(r) = max(maxval(abs(mominv(1:21))),maxval(abs(masses2(0:6))))
5074  if(norm(r).ne.0d0) then
5075  norm(r)=1d0/norm(r)**(5-real(r)/2)
5076  else
5077  norm(r)=1d0/muir2_cll**(5-real(r)/2)
5078  end if
5079  end if
5080  tgacc(r) = tgerr_aux(r)/norm(r)
5081  end do
5082 
5083  end if
5084 
5085  else
5086  call g_main_cll(cg,cguv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5),mominv(6), &
5087  mominv(7),mominv(8),mominv(9),mominv(10),mominv(11),mominv(12), &
5088  mominv(13),mominv(14),mominv(15),mominv(16),mominv(17),mominv(18), &
5089  mominv(19),mominv(20),mominv(21),masses2(0),masses2(1),masses2(2), &
5090  masses2(3),masses2(4),masses2(5),masses2(6),rmax,gerr2=cgerr,id_in=0)
5091  call calctensorg(tg,tguv,tgerr_aux,cg,cguv,cgerr,momvec,rmax)
5092  if (present(tgerr)) tgerr = tgerr_aux
5093  norm = 0d0
5094  do r=0,rmax
5095  do n0=0,r
5096  do n1=0,r-n0
5097  do n2=0,r-n0-n1
5098  n3=r-n0-n1-n2
5099  norm(r) = max(norm(r),abs(tg(n0,n1,n2,n3)))
5100  end do
5101  end do
5102  end do
5103  if (norm(r).eq.0d0) then
5104  norm(r) = max(maxval(abs(mominv(1:21))),maxval(abs(masses2(0:6))))
5105  if(norm(r).ne.0d0) then
5106  norm(r)=1d0/norm(r)**(5-real(r)/2)
5107  else
5108  norm(r)=1d0/muir2_cll**(5-real(r)/2)
5109  end if
5110  end if
5111  tgacc(r) = tgerr_aux(r)/norm(r)
5112  end do
5113  end if
5114 
5115  if (monitoring) then
5116  if (monitoring) pointscntgten_cll = pointscntgten_cll + 1
5117 
5118  if(maxval(tgacc).gt.reqacc_cll) accpointscntgten_cll = accpointscntgten_cll + 1
5119 
5120  if(maxval(tgacc).gt.critacc_cll) then
5121  critpointscntgten_cll = critpointscntgten_cll + 1
5122  if ( critpointscntgten_cll.le.noutcritpointsmax_cll(7) ) then
5123  call critpointsout_cll('TGten_cll',0,maxval(tgacc),critpointscntgten_cll)
5124  if( critpointscntgten_cll.eq.noutcritpointsmax_cll(7)) then
5125  write(ncpout_cll,*) ' Further output of Critical Points for TGten_cll suppressed'
5126  write(ncpout_cll,*)
5127  endif
5128  end if
5129  end if
5130  end if
5131 
5132  if (monitoring) then
5133  pointscntgten_cll = pointscntgten_cll + 1
5134 
5135  if(maxval(tgacc).gt.reqacc_cll) accpointscntgten_cll = accpointscntgten_cll + 1
5136 
5137  if(maxval(tgacc).gt.critacc_cll) then
5138  critpointscntgten_cll = critpointscntgten_cll + 1
5139  if ( critpointscntgten_cll.le.noutcritpointsmax_cll(7) ) then
5140  call critpointsout_cll('TGten_cll',0,maxval(tgacc),critpointscntgten_cll)
5141  if( critpointscntgten_cll.eq.noutcritpointsmax_cll(7)) then
5142  write(ncpout_cll,*) ' Further output of Critical Points for TGten_cll suppressed'
5143  write(ncpout_cll,*)
5144  endif
5145 #ifdef CritPoints2
5146  call critpointsout2_cll('TGten_cll',0,maxval(tgacc),critpointscntgten_cll)
5147  if( critpointscntgten_cll.eq.noutcritpointsmax_cll(7)) then
5148  write(ncpout2_cll,*) ' Further output of Critical Points for TGten_cll suppressed'
5149  write(ncpout2_cll,*)
5150  endif
5151 #endif
5152  end if
5153  end if
5154  end if
5155 

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