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

Public Member Functions

subroutine ften_main_cll (TF, TFuv, MomVec, MomInv, masses2, rmax, TFerr)
 
subroutine ften_list_cll (TF, TFuv, MomVec, MomInv, masses2, rmax, TFerr)
 
subroutine ften_args_cll (TF, TFuv, p1vec, p2vec, p3vec, p4vec, p5vec, p10, p21, p32, p43, p54, p50, p20, p31, p42, p53, p40, p51, p30, p41, p52, m02, m12, m22, m32, m42, m52, rmax, TFerr)
 
subroutine ften_args_list_cll (TF, TFuv, p1vec, p2vec, p3vec, p4vec, p5vec, p10, p21, p32, p43, p54, p50, p20, p31, p42, p53, p40, p51, p30, p41, p52, m02, m12, m22, m32, m42, m52, rmax, TFerr)
 

Detailed Description

Definition at line 74 of file collier_tensors.F90.

Member Function/Subroutine Documentation

◆ ften_args_cll()

subroutine collier_tensors::ften_cll::ften_args_cll ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TF,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TFuv,
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, 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)  p50,
double complex, intent(in)  p20,
double complex, intent(in)  p31,
double complex, intent(in)  p42,
double complex, intent(in)  p53,
double complex, intent(in)  p40,
double complex, intent(in)  p51,
double complex, intent(in)  p30,
double complex, intent(in)  p41,
double complex, intent(in)  p52,
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,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  TFerr 
)

Definition at line 4325 of file collier_tensors.F90.

4325  integer, intent(in) :: rmax
4326  double complex, intent(in) :: p1vec(0:3),p2vec(0:3),p3vec(0:3),p4vec(0:3)
4327  double complex, intent(in) :: p5vec(0:3)
4328  double complex, intent(in) :: p10,p21,p32,p43,p54,p50,p20,p31,p42,p53,p40
4329  double complex, intent(in) :: p51,p30,p41,p52,m02,m12,m22,m32,m42,m52
4330  double complex, intent(out) :: TF(0:rmax,0:rmax,0:rmax,0:rmax)
4331  double complex, intent(out) :: TFuv(0:rmax,0:rmax,0:rmax,0:rmax)
4332  double precision, intent(out), optional :: TFerr(0:rmax)
4333  double complex :: TF2(0:rmax,0:rmax,0:rmax,0:rmax), TFuv2(0:rmax,0:rmax,0:rmax,0:rmax)
4334  double complex :: MomVec(0:3,5), MomInv(15), masses2(0:5)
4335  double complex :: CF(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
4336  double complex :: CFuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
4337  double precision :: CFerr(0:rmax), TFerr_aux(0:rmax), TFerr_aux2(0:rmax)
4338  double complex :: args(41)
4339  double precision :: TFdiff(0:rmax),norm(0:rmax),norm_coli,norm_dd,TFacc(0:rmax)
4340  integer :: r,n0,n1,n2,n3
4341  logical :: eflag
4342 
4343  if (6.gt.nmax_cll) then
4344  call seterrflag_cll(-10)
4345  call errout_cll('Ften_cll','Nmax_cll smaller 6',eflag,.true.)
4346  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
4347  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 6'
4348  call propagateerrflag_cll
4349  return
4350  end if
4351  if (rmax.gt.rmax_cll) then
4352  call seterrflag_cll(-10)
4353  call errout_cll('Ften_cll','argument rmax larger than rmax_cll',eflag,.true.)
4354  write(nerrout_cll,*) 'rmax =',rmax
4355  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
4356  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
4357  call propagateerrflag_cll
4358  return
4359  end if
4360 
4361  momvec(0:,1) = p1vec
4362  momvec(0:,2) = p2vec
4363  momvec(0:,3) = p3vec
4364  momvec(0:,4) = p4vec
4365  momvec(0:,5) = p5vec
4366  mominv(1) = p10
4367  mominv(2) = p21
4368  mominv(3) = p32
4369  mominv(4) = p43
4370  mominv(5) = p54
4371  mominv(6) = p50
4372  mominv(7) = p20
4373  mominv(8) = p31
4374  mominv(9) = p42
4375  mominv(10) = p53
4376  mominv(11) = p40
4377  mominv(12) = p51
4378  mominv(13) = p30
4379  mominv(14) = p41
4380  mominv(15) = p52
4381  masses2(0) = m02
4382  masses2(1) = m12
4383  masses2(2) = m22
4384  masses2(3) = m32
4385  masses2(4) = m42
4386  masses2(5) = m52
4387 
4388  ! set ID of master call
4389  args(1:4) = momvec(0:,1)
4390  args(5:8) = momvec(0:,2)
4391  args(9:12) = momvec(0:,3)
4392  args(13:16) = momvec(0:,4)
4393  args(17:20) = momvec(0:,5)
4394  args(21:35) = mominv
4395  args(36:41) = masses2(0:)
4396  call setmasterfname_cll('Ften_cll')
4397  call setmastern_cll(6)
4398  call setmasterr_cll(rmax)
4399  call setmasterargs_cll(41,args)
4400 
4401  call settencache_cll(tenred_cll-1)
4402 
4403  if (tenred_cll.le.6) then
4404 
4405  if (mode_cll.gt.1) call f_dd_dummy(rmax)
4406 
4407  if (mode_cll.eq.3) then
4408  ! calculate tensor with coefficients from COLI
4409  mode_cll = 1
4410  call calctensorfr(tf,tfuv,tferr_aux,momvec,mominv,masses2,rmax)
4411 
4412  ! calculate tensor with coefficients from DD
4413  mode_cll = 2
4414  call calctensorfr(tf2,tfuv2,tferr_aux2,momvec,mominv,masses2,rmax)
4415 
4416  ! comparison --> take better result
4417  mode_cll = 3
4418  do r=0,rmax
4419  norm_coli=0d0
4420  norm_dd=0d0
4421  do n0=0,r
4422  do n1=0,r-n0
4423  do n2=0,r-n0-n1
4424  n3=r-n0-n1-n2
4425  norm_coli = max(norm_coli,abs(tf(n0,n1,n2,n3)))
4426  norm_dd = max(norm_dd,abs(tf2(n0,n1,n2,n3)))
4427  end do
4428  end do
4429  end do
4430  if (norm_coli.eq.0d0) then
4431  norm_coli = max(maxval(abs(mominv(1:15))),maxval(abs(masses2(0:5))))
4432  if(norm_coli.ne.0d0) then
4433  norm_coli=1d0/norm_coli**(4-real(r)/2)
4434  else
4435  norm_coli=1d0/muir2_cll**(4-real(r)/2)
4436  end if
4437  end if
4438  if (norm_dd.eq.0d0) then
4439  norm_dd = max(maxval(abs(mominv(1:15))),maxval(abs(masses2(0:5))))
4440  if(norm_dd.ne.0d0) then
4441  norm_dd=1d0/norm_dd**(4-real(r)/2)
4442  else
4443  norm_dd=1d0/muir2_cll**(4-real(r)/2)
4444  end if
4445  end if
4446  norm(r) = min(norm_coli,norm_dd)
4447  end do
4448 
4449  call checktensors_cll(tf,tf2,momvec,mominv,masses2,norm,6,rmax,tfdiff)
4450 
4451  if (tferr_aux(rmax).lt.tferr_aux2(rmax)) then
4452  if (present(tferr)) tferr = max(tferr_aux,tfdiff*norm)
4453  do r=0,rmax
4454  tfacc(r) = max(tferr_aux(r)/norm(r),tfdiff(r))
4455  end do
4456  if (monitoring) pointscntften_coli = pointscntften_coli + 1
4457  else
4458  tf = tf2
4459  tfuv = tfuv2
4460  if (present(tferr)) tferr = max(tferr_aux2,tfdiff*norm)
4461  do r=0,rmax
4462  tfacc(r) = max(tferr_aux2(r)/norm(r),tfdiff(r))
4463  end do
4464  if (monitoring) pointscntften_dd = pointscntften_dd + 1
4465  end if
4466 
4467  else
4468  call calctensorfr(tf,tfuv,tferr_aux,momvec,mominv,masses2,rmax)
4469  if (present(tferr)) tferr = tferr_aux
4470  norm = 0d0
4471  do r=0,rmax
4472  do n0=0,r
4473  do n1=0,r-n0
4474  do n2=0,r-n0-n1
4475  n3=r-n0-n1-n2
4476  norm(r) = max(norm(r),abs(tf(n0,n1,n2,n3)))
4477  end do
4478  end do
4479  end do
4480  if (norm(r).eq.0d0) then
4481  norm(r) = max(maxval(abs(mominv(1:15))),maxval(abs(masses2(0:5))))
4482  if(norm(r).ne.0d0) then
4483  norm(r)=1d0/norm(r)**(4-real(r)/2)
4484  else
4485  norm(r)=1d0/muir2_cll**(4-real(r)/2)
4486  end if
4487  end if
4488  tfacc(r) = tferr_aux(r)/norm(r)
4489  end do
4490 
4491  end if
4492 
4493  else
4494 
4495  if (mode_cll.eq.3) then
4496  ! calculate tensor with coefficients from COLI
4497  mode_cll = 1
4498  call f_main_cll(cf,cfuv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5), &
4499  mominv(6),mominv(7),mominv(8),mominv(9),mominv(10),mominv(11), &
4500  mominv(12),mominv(13),mominv(14),mominv(15),masses2(0),masses2(1), &
4501  masses2(2),masses2(3),masses2(4),masses2(5),rmax,ferr2=cferr,id_in=0)
4502  call calctensorf(tf,tfuv,tferr_aux,cf,cfuv,cferr,momvec,rmax)
4503 
4504  ! calculate tensor with coefficients from DD
4505  mode_cll = 2
4506  call f_main_cll(cf,cfuv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5), &
4507  mominv(6),mominv(7),mominv(8),mominv(9),mominv(10),mominv(11), &
4508  mominv(12),mominv(13),mominv(14),mominv(15),masses2(0),masses2(1), &
4509  masses2(2),masses2(3),masses2(4),masses2(5),rmax,ferr2=cferr,id_in=0)
4510  call calctensorf(tf2,tfuv2,tferr_aux2,cf,cfuv,cferr,momvec,rmax)
4511 
4512  ! comparison --> take better result
4513  mode_cll = 3
4514  do r=0,rmax
4515  norm_coli=0d0
4516  norm_dd=0d0
4517  do n0=0,r
4518  do n1=0,r-n0
4519  do n2=0,r-n0-n1
4520  n3=r-n0-n1-n2
4521  norm_coli = max(norm_coli,abs(tf(n0,n1,n2,n3)))
4522  norm_dd = max(norm_dd,abs(tf2(n0,n1,n2,n3)))
4523  end do
4524  end do
4525  end do
4526  if (norm_coli.eq.0d0) then
4527  norm_coli = max(maxval(abs(mominv(1:15))),maxval(abs(masses2(0:5))))
4528  if(norm_coli.ne.0d0) then
4529  norm_coli=1d0/norm_coli**(4-real(r)/2)
4530  else
4531  norm_coli=1d0/muir2_cll**(4-real(r)/2)
4532  end if
4533  end if
4534  if (norm_dd.eq.0d0) then
4535  norm_dd = max(maxval(abs(mominv(1:15))),maxval(abs(masses2(0:5))))
4536  if(norm_dd.ne.0d0) then
4537  norm_dd=1d0/norm_dd**(4-real(r)/2)
4538  else
4539  norm_dd=1d0/muir2_cll**(4-real(r)/2)
4540  end if
4541  end if
4542  norm(r) = min(norm_coli,norm_dd)
4543  end do
4544 
4545  call checktensors_cll(tf,tf2,momvec,mominv,masses2,norm,6,rmax,tfdiff)
4546 
4547  if (tferr_aux(rmax).lt.tferr_aux2(rmax)) then
4548  if (present(tferr)) tferr = max(tferr_aux,tfdiff*norm)
4549  do r=0,rmax
4550  tfacc(r) = max(tferr_aux(r)/norm(r),tfdiff(r))
4551  end do
4552  if (monitoring) pointscntften_coli = pointscntften_coli + 1
4553  else
4554  tf = tf2
4555  tfuv = tfuv2
4556  if (present(tferr)) tferr = max(tferr_aux2,tfdiff*norm)
4557  do r=0,rmax
4558  tfacc(r) = max(tferr_aux2(r)/norm(r),tfdiff(r))
4559  end do
4560  if (monitoring) pointscntften_dd = pointscntften_dd + 1
4561  end if
4562 
4563  else
4564  call f_main_cll(cf,cfuv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5), &
4565  mominv(6),mominv(7),mominv(8),mominv(9),mominv(10),mominv(11), &
4566  mominv(12),mominv(13),mominv(14),mominv(15),masses2(0),masses2(1), &
4567  masses2(2),masses2(3),masses2(4),masses2(5),rmax,ferr2=cferr,id_in=0)
4568  call calctensorf(tf,tfuv,tferr_aux,cf,cfuv,cferr,momvec,rmax)
4569  if (present(tferr)) tferr = tferr_aux
4570  norm = 0d0
4571  do r=0,rmax
4572  do n0=0,r
4573  do n1=0,r-n0
4574  do n2=0,r-n0-n1
4575  n3=r-n0-n1-n2
4576  norm(r) = max(norm(r),abs(tf(n0,n1,n2,n3)))
4577  end do
4578  end do
4579  end do
4580  if (norm(r).eq.0d0) then
4581  norm(r) = max(maxval(abs(mominv(1:15))),maxval(abs(masses2(0:5))))
4582  if(norm(r).ne.0d0) then
4583  norm(r)=1d0/norm(r)**(4-real(r)/2)
4584  else
4585  norm(r)=1d0/muir2_cll**(4-real(r)/2)
4586  end if
4587  end if
4588  tfacc(r) = tferr_aux(r)/norm(r)
4589  end do
4590 
4591  end if
4592 
4593  end if
4594 
4595  call propagateaccflag_cll(tfacc,rmax)
4596  call propagateerrflag_cll
4597 
4598  if (monitoring) then
4599  pointscntften_cll = pointscntften_cll + 1
4600 
4601  if(maxval(tfacc).gt.reqacc_cll) accpointscntften_cll = accpointscntften_cll + 1
4602 
4603  if(maxval(tfacc).gt.critacc_cll) then
4604  critpointscntften_cll = critpointscntften_cll + 1
4605  if ( critpointscntften_cll.le.noutcritpointsmax_cll(6) ) then
4606  call critpointsout_cll('TFten_cll',0,maxval(tfacc),critpointscntften_cll)
4607  if( critpointscntften_cll.eq.noutcritpointsmax_cll(6)) then
4608  write(ncpout_cll,*) ' Further output of Critical Points for TFten_cll suppressed'
4609  write(ncpout_cll,*)
4610  endif
4611 #ifdef CritPoints2
4612  call critpointsout2_cll('TFten_cll',0,maxval(tfacc),critpointscntften_cll)
4613  if( critpointscntften_cll.eq.noutcritpointsmax_cll(6)) then
4614  write(ncpout2_cll,*) ' Further output of Critical Points for TFten_cll suppressed'
4615  write(ncpout2_cll,*)
4616  endif
4617 #endif
4618  end if
4619  end if
4620  end if
4621 

◆ ften_args_list_cll()

subroutine collier_tensors::ften_cll::ften_args_list_cll ( double complex, dimension(:), intent(out)  TF,
double complex, dimension(:), intent(out)  TFuv,
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, 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)  p50,
double complex, intent(in)  p20,
double complex, intent(in)  p31,
double complex, intent(in)  p42,
double complex, intent(in)  p53,
double complex, intent(in)  p40,
double complex, intent(in)  p51,
double complex, intent(in)  p30,
double complex, intent(in)  p41,
double complex, intent(in)  p52,
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,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  TFerr 
)

Definition at line 4638 of file collier_tensors.F90.

4638  integer, intent(in) :: rmax
4639  double complex, intent(in) :: p1vec(0:3),p2vec(0:3),p3vec(0:3),p4vec(0:3)
4640  double complex, intent(in) :: p5vec(0:3)
4641  double complex, intent(in) :: p10,p21,p32,p43,p54,p50,p20,p31,p42,p53,p40
4642  double complex, intent(in) :: p51,p30,p41,p52,m02,m12,m22,m32,m42,m52
4643  double complex, intent(out) :: TF(:),TFuv(:)
4644  double precision, intent(out), optional :: TFerr(0:rmax)
4645  logical :: eflag
4646 
4647  if (6.gt.nmax_cll) then
4648  call seterrflag_cll(-10)
4649  call errout_cll('Ften_cll','Nmax_cll smaller 6',eflag,.true.)
4650  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
4651  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 6'
4652  call propagateerrflag_cll
4653  return
4654  end if
4655  if (rmax.gt.rmax_cll) then
4656  call seterrflag_cll(-10)
4657  call errout_cll('Ften_cll','argument rmax larger than rmax_cll',eflag,.true.)
4658  write(nerrout_cll,*) 'rmax =',rmax
4659  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
4660  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
4661  call propagateerrflag_cll
4662  return
4663  end if
4664 
4665  call ften_args_list_checked_cll(tf,tfuv,p1vec,p2vec,p3vec,p4vec,p5vec, &
4666  p10,p21,p32,p43,p54,p50,p20,p31,p42,p53,p40, &
4667  p51,p30,p41,p52,m02,m12,m22,m32,m42,m52,rmax,tferr)
4668 

◆ ften_list_cll()

subroutine collier_tensors::ften_cll::ften_list_cll ( double complex, dimension(:), intent(out)  TF,
double complex, dimension(:), intent(out)  TFuv,
double complex, dimension(0:3,5), intent(in)  MomVec,
double complex, dimension(15), intent(in)  MomInv,
double complex, dimension(0:5), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  TFerr 
)

Definition at line 4049 of file collier_tensors.F90.

4049 
4050  integer, intent(in) :: rmax
4051  double complex, intent(in) :: MomVec(0:3,5), MomInv(15), masses2(0:5)
4052  double complex, intent(out) :: TF(:), TFuv(:)
4053  double precision, intent(out), optional :: TFerr(0:rmax)
4054  logical :: eflag
4055 
4056  if (6.gt.nmax_cll) then
4057  call seterrflag_cll(-10)
4058  call errout_cll('Ften_cll','Nmax_cll smaller 6',eflag,.true.)
4059  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
4060  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 6'
4061  call propagateerrflag_cll
4062  return
4063  end if
4064  if (rmax.gt.rmax_cll) then
4065  call seterrflag_cll(-10)
4066  call errout_cll('Ften_cll','argument rmax larger than rmax_cll',eflag,.true.)
4067  write(nerrout_cll,*) 'rmax =',rmax
4068  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
4069  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
4070  call propagateerrflag_cll
4071  return
4072  end if
4073 
4074  call ften_list_checked_cll(tf,tfuv,momvec,mominv,masses2,rmax,tferr)
4075 

◆ ften_main_cll()

subroutine collier_tensors::ften_cll::ften_main_cll ( double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TF,
double complex, dimension(0:rmax,0:rmax,0:rmax,0:rmax), intent(out)  TFuv,
double complex, dimension(0:3,5), intent(in)  MomVec,
double complex, dimension(15), intent(in)  MomInv,
double complex, dimension(0:5), intent(in)  masses2,
integer, intent(in)  rmax,
double precision, dimension(0:rmax), intent(out), optional  TFerr 
)

Definition at line 3768 of file collier_tensors.F90.

3768 
3769  integer, intent(in) :: rmax
3770  double complex, intent(in) :: MomVec(0:3,5), MomInv(15), masses2(0:5)
3771  double complex, intent(out) :: TF(0:rmax,0:rmax,0:rmax,0:rmax)
3772  double complex, intent(out) :: TFuv(0:rmax,0:rmax,0:rmax,0:rmax)
3773  double precision, intent(out), optional :: TFerr(0:rmax)
3774  double complex :: TF2(0:rmax,0:rmax,0:rmax,0:rmax), TFuv2(0:rmax,0:rmax,0:rmax,0:rmax)
3775  double complex :: CF(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
3776  double complex :: CFuv(0:rmax/2,0:rmax,0:rmax,0:rmax,0:rmax,0:rmax)
3777  double precision :: CFerr(0:rmax), TFerr_aux(0:rmax), TFerr_aux2(0:rmax)
3778  double complex :: args(41)
3779  double precision :: TFdiff(0:rmax),norm(0:rmax),norm_coli,norm_dd,TFacc(0:rmax)
3780  integer :: r,n0,n1,n2,n3
3781  logical :: eflag
3782 
3783  if (6.gt.nmax_cll) then
3784  call seterrflag_cll(-10)
3785  call errout_cll('Ften_cll','Nmax_cll smaller 6',eflag,.true.)
3786  write(nerrout_cll,*) 'Nmax_cll =',nmax_cll
3787  write(nerrout_cll,*) 'Reinitialize COLLIER with Nmax_cll >= 6'
3788  call propagateerrflag_cll
3789  return
3790  end if
3791  if (rmax.gt.rmax_cll) then
3792  call seterrflag_cll(-10)
3793  call errout_cll('Ften_cll','argument rmax larger than rmax_cll',eflag,.true.)
3794  write(nerrout_cll,*) 'rmax =',rmax
3795  write(nerrout_cll,*) 'rmax_cll =',rmax_cll
3796  write(nerrout_cll,*) 'Reinitialize COLLIER with rmax_cll >= ',rmax
3797  call propagateerrflag_cll
3798  return
3799  end if
3800 
3801  ! set ID of master call
3802  args(1:4) = momvec(0:,1)
3803  args(5:8) = momvec(0:,2)
3804  args(9:12) = momvec(0:,3)
3805  args(13:16) = momvec(0:,4)
3806  args(17:20) = momvec(0:,5)
3807  args(21:35) = mominv
3808  args(36:41) = masses2(0:)
3809  call setmasterfname_cll('Ften_cll')
3810  call setmastern_cll(6)
3811  call setmasterr_cll(rmax)
3812  call setmasterargs_cll(41,args)
3813 
3814  call settencache_cll(tenred_cll-1)
3815 
3816 
3817  if (tenred_cll.le.6) then
3818 
3819  if (mode_cll.gt.1) call f_dd_dummy(rmax)
3820 
3821  if (mode_cll.eq.3) then
3822  ! calculate tensor with coefficients from COLI
3823  mode_cll = 1
3824  call calctensorfr(tf,tfuv,tferr_aux,momvec,mominv,masses2,rmax)
3825 
3826  ! calculate tensor with coefficients from DD
3827  mode_cll = 2
3828  call calctensorfr(tf2,tfuv2,tferr_aux2,momvec,mominv,masses2,rmax)
3829 
3830  ! comparison --> take better result
3831  mode_cll = 3
3832  do r=0,rmax
3833  norm_coli=0d0
3834  norm_dd=0d0
3835  do n0=0,r
3836  do n1=0,r-n0
3837  do n2=0,r-n0-n1
3838  n3=r-n0-n1-n2
3839  norm_coli = max(norm_coli,abs(tf(n0,n1,n2,n3)))
3840  norm_dd = max(norm_dd,abs(tf2(n0,n1,n2,n3)))
3841  end do
3842  end do
3843  end do
3844  if (norm_coli.eq.0d0) then
3845  norm_coli = max(maxval(abs(mominv(1:15))),maxval(abs(masses2(0:5))))
3846  if(norm_coli.ne.0d0) then
3847  norm_coli=1d0/norm_coli**(4-real(r)/2)
3848  else
3849  norm_coli=1d0/muir2_cll**(4-real(r)/2)
3850  end if
3851  end if
3852  if (norm_dd.eq.0d0) then
3853  norm_dd = max(maxval(abs(mominv(1:15))),maxval(abs(masses2(0:5))))
3854  if(norm_dd.ne.0d0) then
3855  norm_dd=1d0/norm_dd**(4-real(r)/2)
3856  else
3857  norm_dd=1d0/muir2_cll**(4-real(r)/2)
3858  end if
3859  end if
3860  norm(r) = min(norm_coli,norm_dd)
3861  end do
3862 
3863  call checktensors_cll(tf,tf2,momvec,mominv,masses2,norm,6,rmax,tfdiff)
3864 
3865  if (tferr_aux(rmax).lt.tferr_aux2(rmax)) then
3866  if (present(tferr)) tferr = max(tferr_aux,tfdiff*norm)
3867  do r=0,rmax
3868  tfacc(r) = max(tferr_aux(r)/norm(r),tfdiff(r))
3869  end do
3870  if (monitoring) pointscntften_coli = pointscntften_coli + 1
3871  else
3872  tf = tf2
3873  tfuv = tfuv2
3874  if (present(tferr)) tferr = max(tferr_aux2,tfdiff*norm)
3875  do r=0,rmax
3876  tfacc(r) = max(tferr_aux2(r)/norm(r),tfdiff(r))
3877  end do
3878  if (monitoring) pointscntften_dd = pointscntften_dd + 1
3879  end if
3880 
3881  else
3882  call calctensorfr(tf,tfuv,tferr_aux,momvec,mominv,masses2,rmax)
3883  if (present(tferr)) tferr = tferr_aux
3884  norm = 0d0
3885  do r=0,rmax
3886  do n0=0,r
3887  do n1=0,r-n0
3888  do n2=0,r-n0-n1
3889  n3=r-n0-n1-n2
3890  norm(r) = max(norm(r),abs(tf(n0,n1,n2,n3)))
3891  end do
3892  end do
3893  end do
3894  if (norm(r).eq.0d0) then
3895  norm(r) = max(maxval(abs(mominv(1:15))),maxval(abs(masses2(0:5))))
3896  if(norm(r).ne.0d0) then
3897  norm(r)=1d0/norm(r)**(4-real(r)/2)
3898  else
3899  norm(r)=1d0/muir2_cll**(4-real(r)/2)
3900  end if
3901  end if
3902  tfacc(r) = tferr_aux(r)/norm(r)
3903  end do
3904 
3905  end if
3906 
3907 
3908  else
3909 
3910  if (mode_cll.eq.3) then
3911  ! calculate tensor with coefficients from COLI
3912  mode_cll = 1
3913  call f_main_cll(cf,cfuv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5), &
3914  mominv(6),mominv(7),mominv(8),mominv(9),mominv(10),mominv(11), &
3915  mominv(12),mominv(13),mominv(14),mominv(15),masses2(0),masses2(1), &
3916  masses2(2),masses2(3),masses2(4),masses2(5),rmax,ferr2=cferr,id_in=0)
3917  call calctensorf(tf,tfuv,tferr_aux,cf,cfuv,cferr,momvec,rmax)
3918 
3919  ! calculate tensor with coefficients from DD
3920  mode_cll = 2
3921  call f_main_cll(cf,cfuv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5), &
3922  mominv(6),mominv(7),mominv(8),mominv(9),mominv(10),mominv(11), &
3923  mominv(12),mominv(13),mominv(14),mominv(15),masses2(0),masses2(1), &
3924  masses2(2),masses2(3),masses2(4),masses2(5),rmax,ferr2=cferr,id_in=0)
3925  call calctensorf(tf2,tfuv2,tferr_aux2,cf,cfuv,cferr,momvec,rmax)
3926 
3927  ! comparison --> take better result
3928  mode_cll = 3
3929  do r=0,rmax
3930  norm_coli=0d0
3931  norm_dd=0d0
3932  do n0=0,r
3933  do n1=0,r-n0
3934  do n2=0,r-n0-n1
3935  n3=r-n0-n1-n2
3936  norm_coli = max(norm_coli,abs(tf(n0,n1,n2,n3)))
3937  norm_dd = max(norm_dd,abs(tf2(n0,n1,n2,n3)))
3938  end do
3939  end do
3940  end do
3941  if (norm_coli.eq.0d0) then
3942  norm_coli = max(maxval(abs(mominv(1:15))),maxval(abs(masses2(0:5))))
3943  if(norm_coli.ne.0d0) then
3944  norm_coli=1d0/norm_coli**(4-real(r)/2)
3945  else
3946  norm_coli=1d0/muir2_cll**(4-real(r)/2)
3947  end if
3948  end if
3949  if (norm_dd.eq.0d0) then
3950  norm_dd = max(maxval(abs(mominv(1:15))),maxval(abs(masses2(0:5))))
3951  if(norm_dd.ne.0d0) then
3952  norm_dd=1d0/norm_dd**(4-real(r)/2)
3953  else
3954  norm_dd=1d0/muir2_cll**(4-real(r)/2)
3955  end if
3956  end if
3957  norm(r) = min(norm_coli,norm_dd)
3958  end do
3959 
3960  call checktensors_cll(tf,tf2,momvec,mominv,masses2,norm,6,rmax,tfdiff)
3961 
3962  if (tferr_aux(rmax).lt.tferr_aux2(rmax)) then
3963  if (present(tferr)) tferr = max(tferr_aux,tfdiff*norm)
3964  do r=0,rmax
3965  tfacc(r) = max(tferr_aux(r)/norm(r),tfdiff(r))
3966  end do
3967  if (monitoring) pointscntften_coli = pointscntften_coli + 1
3968  else
3969  tf = tf2
3970  tfuv = tfuv2
3971  if (present(tferr)) tferr = max(tferr_aux2,tfdiff*norm)
3972  do r=0,rmax
3973  tfacc(r) = max(tferr_aux2(r)/norm(r),tfdiff(r))
3974  end do
3975  if (monitoring) pointscntften_dd = pointscntften_dd + 1
3976  end if
3977 
3978  else
3979  call f_main_cll(cf,cfuv,mominv(1),mominv(2),mominv(3),mominv(4),mominv(5), &
3980  mominv(6),mominv(7),mominv(8),mominv(9),mominv(10),mominv(11), &
3981  mominv(12),mominv(13),mominv(14),mominv(15),masses2(0),masses2(1), &
3982  masses2(2),masses2(3),masses2(4),masses2(5),rmax,ferr2=cferr,id_in=0)
3983  call calctensorf(tf,tfuv,tferr_aux,cf,cfuv,cferr,momvec,rmax)
3984  if (present(tferr)) tferr = tferr_aux
3985  norm = 0d0
3986  do r=0,rmax
3987  do n0=0,r
3988  do n1=0,r-n0
3989  do n2=0,r-n0-n1
3990  n3=r-n0-n1-n2
3991  norm(r) = max(norm(r),abs(tf(n0,n1,n2,n3)))
3992  end do
3993  end do
3994  end do
3995  if (norm(r).eq.0d0) then
3996  norm(r) = max(maxval(abs(mominv(1:15))),maxval(abs(masses2(0:5))))
3997  if(norm(r).ne.0d0) then
3998  norm(r)=1d0/norm(r)**(4-real(r)/2)
3999  else
4000  norm(r)=1d0/muir2_cll**(4-real(r)/2)
4001  end if
4002  end if
4003  tfacc(r) = tferr_aux(r)/norm(r)
4004  end do
4005 
4006  end if
4007 
4008  end if
4009 
4010  call propagateaccflag_cll(tfacc,rmax)
4011  call propagateerrflag_cll
4012 
4013  if (monitoring) then
4014  pointscntften_cll = pointscntften_cll + 1
4015 
4016  if(maxval(tfacc).gt.reqacc_cll) accpointscntften_cll = accpointscntften_cll + 1
4017 
4018  if(maxval(tfacc).gt.critacc_cll) then
4019  critpointscntften_cll = critpointscntften_cll + 1
4020  if ( critpointscntften_cll.le.noutcritpointsmax_cll(6) ) then
4021  call critpointsout_cll('TFten_cll',0,maxval(tfacc),critpointscntften_cll)
4022  if( critpointscntften_cll.eq.noutcritpointsmax_cll(6)) then
4023  write(ncpout_cll,*) ' Further output of Critical Points for TFten_cll suppressed'
4024  write(ncpout_cll,*)
4025  endif
4026 #ifdef CritPoints2
4027  call critpointsout2_cll('TFten_cll',0,maxval(tfacc),critpointscntften_cll)
4028  if( critpointscntften_cll.eq.noutcritpointsmax_cll(6)) then
4029  write(ncpout2_cll,*) ' Further output of Critical Points for TFten_cll suppressed'
4030  write(ncpout2_cll,*)
4031  endif
4032 #endif
4033  end if
4034  end if
4035  end if
4036 

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