JHUGen MELA  JHUGen v7.5.6, MELA v2.4.2
Matrix element calculations as used in JHUGen.
InitTensors.F90
Go to the documentation of this file.
1 !!
2 !! File InitTensors.F90 is part of COLLIER
3 !! - A Complex One-Loop Library In Extended Regularizations
4 !!
5 !! Copyright (C) 2015, 2016 Ansgar Denner, Stefan Dittmaier, Lars Hofer
6 !!
7 !! COLLIER is licenced under the GNU GPL version 3, see COPYING for details.
8 !!
9 
10 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
11 !
12 ! *************************
13 ! * module InitTensors *
14 ! * by Lars Hofer *
15 ! *************************
16 !
17 ! global variables:
18 !
19 ! functions and subroutines:
20 !
21 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
22 
23 
24 
26 
27  use combinatorics
28 
29  implicit none
30 
31  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
32  ! global variable CFtab
33  !
34  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
35 
36 
37  integer, dimension(:,:), allocatable :: cftab, lorindtab, addindtab, addgtab
38  integer, dimension(:,:,:,:), allocatable :: indspiprod
39  integer, dimension(:), allocatable :: rts
40  real :: stin, stout, time1, time2, time3
41  logical :: calcuv_ten
42 
43 
44 contains
45 
46  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
47  ! subroutine init_tables2(rmax)
48  !
49  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
50 
51  subroutine init_tables2(Nm1,rmax)
52 
53  integer, intent(in) :: Nm1,rmax
54 
55  call setrts(rmax)
56  call setlorindtab(rmax)
57  call setaddindtab(rmax)
58  call setaddgtab(rmax)
59  call setcftab(rmax)
60  call setindspiprod(nm1,rmax)
61 
62  end subroutine init_tables2
63 
64 
65 
66 
67 
68  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
69  ! subroutine SetRtS(rmax)
70  !
71  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
72 
73  subroutine setrts(rmax)
74 
75  integer, intent(in) :: rmax
76  integer :: r
77 
78  if (allocated(rts)) deallocate(rts)
79  allocate(rts(-1:rmax+1))
80 
81  rts(-1) = 0
82  do r=0,rmax+1
83  rts(r) = rts(r-1) + binomtable(r,r+3)
84  end do
85 
86  end subroutine setrts
87 
88 
89 
90 
91 
92  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
93  ! subroutine SetLorIndTab(rmax)
94  !
95  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
96 
97  subroutine setlorindtab(rmax)
98 
99  integer, intent(in) :: rmax
100  integer :: r,mu0,mu1,mu2,mu3,cnt
101 
102  if (allocated(lorindtab)) deallocate(lorindtab)
103  allocate(lorindtab(0:3,rts(rmax)))
104 
105  cnt = 1
106  do r=0,rmax
107  do mu0=r,0,-1
108  do mu1=r-mu0,0,-1
109  do mu2=r-mu0-mu1,0,-1
110  mu3=r-mu0-mu1-mu2
111 
112  lorindtab(0,cnt) = mu0
113  lorindtab(1,cnt) = mu1
114  lorindtab(2,cnt) = mu2
115  lorindtab(3,cnt) = mu3
116 
117  cnt = cnt+1
118 
119  end do
120  end do
121  end do
122  end do
123 
124  end subroutine setlorindtab
125 
126 
127 
128 
129 
130  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
131  ! subroutine SetAddIndTab(rmax)
132  !
133  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
134 
135  subroutine setaddindtab(rmax)
136 
137  integer, intent(in) :: rmax
138  integer :: r,mu,nu,IndMu(0:3),i
139 
140  if (allocated(addindtab)) deallocate(addindtab)
141  allocate(addindtab(rts(rmax-1),0:3))
142 
143  addindtab(1,:) = (/ 2,3,4,5 /)
144  do r=1,rmax-1
145  do mu=rts(r-1)+1,rts(r)
146  indmu = lorindtab(:,mu)
147  do i=0,3
148  indmu(i) = indmu(i)+1
149  do nu=rts(r)+1,rts(r+1)
150  if (all(lorindtab(:,nu).eq.indmu)) then
151  addindtab(mu,i) = nu
152  end if
153  end do
154  indmu(i) = indmu(i)-1
155  end do
156  end do
157  end do
158 
159  end subroutine setaddindtab
160 
161 
162 
163 
164 
165  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
166  ! subroutine SetAddGtab(rmax)
167  !
168  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
169 
170  subroutine setaddgtab(rmax)
171 
172  integer, intent(in) :: rmax
173  integer :: r,nsum,mu,nu,munu,IndMu(0:3),IndN(0:3)
174 
175  if (rmax.le.1) return
176 
177  if (allocated(addgtab)) deallocate(addgtab)
178  allocate(addgtab(rts(rmax),rts(rmax/2)))
179 
180  do nsum=1,rmax/2
181  do nu=rts(nsum-1)+1,rts(nsum)
182  indn = lorindtab(:,nu)
183 
184  do munu=rts(2*nsum-1)+1,rts(2*nsum)
185  if (all(lorindtab(:,munu).eq.2*indn)) then
186  addgtab(1,nu) = munu
187  end if
188  end do
189 
190  do r=1,rmax-2*nsum
191  do mu=rts(r-1)+1,rts(r)
192  indmu = lorindtab(:,mu)
193 
194  do munu=rts(r+2*nsum-1)+1,rts(r+2*nsum)
195 
196  if (all(lorindtab(:,munu).eq.indmu+2*indn)) then
197  addgtab(mu,nu) = munu
198  end if
199 
200  end do
201 
202  end do
203  end do
204 
205  end do
206  end do
207 
208  end subroutine setaddgtab
209 
210 
211 
212 
213 
214  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
215  ! subroutine SetCFtab(rmax)
216  !
217  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
218 
219  subroutine setcftab(rmax)
220 
221  integer, intent(in) :: rmax
222  integer :: mu, nu, IndMu(0:3), IndNu(0:3), CF, i, r
223 
224  if (allocated(cftab)) deallocate(cftab)
225  allocate(cftab(rts(rmax),2:rts(rmax/2)))
226 
227  do r=0,rmax
228  do mu=rts(r-1)+1,rts(r)
229  indmu = lorindtab(:,mu)
230 
231  do nu=2,rts((rmax-r)/2)
232  indnu = lorindtab(:,nu)
233 
234  cf = (-1)**(indnu(1)+indnu(2)+indnu(3))
235  do i=0,3
236  cf = cf*calcfactorial(indmu(i)+2*indnu(i)) &
237  /(2**indnu(i)*calcfactorial(indmu(i))*calcfactorial(indnu(i)))
238  end do
239  cftab(mu,nu) = cf
240 
241  end do
242  end do
243  end do
244 
245  end subroutine setcftab
246 
247 
248 
249 
250 
251  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
252  ! subroutine SetIndsPiProd(Nmax,rmax)
253  !
254  ! sets the global variable IndsPiProd to the table
255  ! IndsPiProd = (CIPP(1,rmax),...,CIPP(Nmax,rmax))
256  ! CIPP = CalcIndsPiProd
257  !
258  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
259 
260  subroutine setindspiprod(Nm1max,rmax)
261 
262  integer, intent(in) :: Nm1max, rmax
263  integer :: strmax, Nm1
264 
265  strmax = binomtable(rmax,rmax+nm1max-1)
266  if (allocated(indspiprod)) then
267  deallocate(indspiprod)
268  end if
269  allocate(indspiprod(0:1,nm1max,strmax,nm1max))
270  indspiprod = 0
271 
272  do nm1=1,nm1max
273  strmax = binomtable(rmax,rmax+nm1-1)
274  indspiprod(0:1,1:nm1,1:strmax,nm1) = calcindspiprod(nm1,rmax)
275  end do
276 
277  end subroutine setindspiprod
278 
279 
280 
281 
282 
283  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
284  ! function CalcIndsPiProd(Nm1,r)
285  !
286  ! calculates resulting indices for product of momentum
287  ! tensor (p_i1,...,p_ik) with p_j (im,j = 1,...,Nm1).
288  ! Further the number of indices equal to j in the final
289  ! result is given.
290  !
291  ! example: Nm1=2, r=3
292  ! momentum tensors:
293  ! rank3: 1: p1p1p1, 2: p1p1p2, 3: p1p2p2, 4: p2p2p2
294  ! rank4: 1: p1p1p1p1, 2: p1p1p1p2, 3: p1p1p2p2,
295  ! 4: p1p2p2p2, 5: p2p2p2p2
296  !
297  ! IndsPiProd(2,1,0:1) = (2,3), IndsPiProd(2,2,0:1) = (3,2)
298  !
299  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
300 
301  function calcindspiprod(Nm1,r) result(IndsPiProd)
302 
303  integer, intent(in) :: nm1, r
304  integer, dimension(:,:,:), allocatable :: indspiprod
305  integer :: str, i, j, inds(r), indsp1(r+1), pos, cnt
306 
307  str = binomtable(r,r+nm1-1)
308  allocate(indspiprod(0:1,1:nm1,1:str))
309 
310  do i=1,str
311  inds = indcombiseq(1:r,i,r,nm1)
312 
313  pos = 1
314  do j=1,nm1
315 
316  cnt = 1
317  do while (pos.le.r)
318  if (inds(pos)>j) then
319  exit
320  else if (inds(pos).eq.j) then
321  cnt = cnt+1
322  end if
323  pos = pos+1
324  end do
325 
326  indsp1(1:pos-1) = inds(1:pos-1)
327  indsp1(pos) = j
328  indsp1(pos+1:r+1) = inds(pos:r)
329 
330  indspiprod(0,j,i) = calcposindcombiseq(nm1,r+1,indsp1)
331  indspiprod(1,j,i) = cnt
332 
333  end do
334 
335  end do
336 
337  end function calcindspiprod
338 
339 
340 
341 
342 
343  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
344  ! subroutine SwitchOffCalcUV_cll()
345  !
346  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
347 
348  subroutine switchoffcalcuv_ten()
349 
350  calcuv_ten = .false.
351 
352  end subroutine switchoffcalcuv_ten
353 
354 
355 
356 
357 
358  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
359  ! subroutine SwitchOnCalcUV_ten()
360  !
361  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
362 
363  subroutine switchoncalcuv_ten()
364 
365  calcuv_ten = .true.
366 
367  end subroutine switchoncalcuv_ten
368 
369 
370 
371 
372 
373  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
374  ! subroutine CheckTensors_cll(T1,T2,MomVec,MomInv,masses2,N,rmax,Tdiff)
375  !
376  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
377 
378  subroutine checktensors_cll(T1,T2,MomVec,MomInv,masses2,norm,N,rmax,Tdiff)
379 
380  integer, intent(in) :: N, rmax
381  double complex, intent(in) :: T1(0:rmax,0:rmax,0:rmax,0:rmax)
382  double complex, intent(in) :: T2(0:rmax,0:rmax,0:rmax,0:rmax)
383  double complex, intent(in) :: MomVec(0:3,N-1)
384  double complex, intent(in) :: MomInv(BinomTable(2,N)), masses2(0:N-1)
385  double precision, intent(in) :: norm(0:rmax)
386  double precision, intent(out) :: Tdiff(0:rmax)
387  integer :: r,n0,n1,n2,n3,i,j,flag
388  double complex :: diffTN
389  double precision :: ratio
390  character(len=*),parameter :: fmt1 = "(A8,'(',i2,') = dcmplx(',d25.18,' , ',d25.18,' )')"
391  character(len=*),parameter :: fmt2 = &
392  "(A6,' TNten(',i1,',',i1,',',i1,',',i1,') = (',d23.16,' , ',d23.16,' ), r=',i2)"
393  character(len=*),parameter :: fmt3 = &
394  "(A7,'(',i1,',',i2,') = dcmplx(',d25.18,' , ',d25.18,' )')"
395 
396 
397  checkcntten_cll(n) = checkcntten_cll(n) + 1
398 
399 ! data CheckCntTN /0/
400 
401  flag=1
402  if(diffcntten_cll(n).ge.maxcheck_cll(n)) flag=0
403  if(ncheckout_cll.eq.-1) flag=0
404 
405  ratio=0d0
406 
407  tdiff=0d0
408  do r=0,rmax
409  do n0=0,r
410  do n1=0,r-n0
411  do n2=0,r-n0-n1
412  n3=r-n0-n1-n2
413 
414  difftn = t1(n0,n1,n2,n3)-t2(n0,n1,n2,n3)
415  tdiff(r) = max(tdiff(r),abs(difftn)/norm(r))
416  if ((abs(difftn).gt.checkacc_cll*norm(r)) .and.(flag.eq.1)) then
417  write(ncheckout_cll,*) '***************************************************************************'
418  write(ncheckout_cll,'(A15,I2,A16,I10)') 'TNten with N =',n,', difference NO.', diffcntten_cll(n)+1
419  write(ncheckout_cll,*) 'COLI and DD do not agree! checkacc =', checkacc_cll
420  write(ncheckout_cll,'(A24,I2,A10,I2,A4,I2)') 'TNten integral with N =', n, ' and rank ', r,' of ',rmax
421  write(ncheckout_cll,*) '---------------------------------------------------------------------------'
422  write(ncheckout_cll,*) 'GLOBAL PARAMETERS:'
423  write(ncheckout_cll,*) 'mode ', mode_cll
424  write(ncheckout_cll,*) 'muUV2 ', muuv2_cll
425  write(ncheckout_cll,*) 'muIR2 ', muir2_cll
426  write(ncheckout_cll,*) 'deltaUV ', deltauv_cll
427  write(ncheckout_cll,*) 'deltaIR1 ', deltair1_cll
428  write(ncheckout_cll,*) 'deltaIR2 ', deltair2_cll
429  write(ncheckout_cll,*) 'nminf ', nminf_cll
430  do i=1,nminf_cll
431  write(ncheckout_cll,*) 'minf2 ', i, minf2_cll(i)
432  end do
433  write(ncheckout_cll,*) 'dprec ', dprec_cll
434  write(ncheckout_cll,*) 'reqacc ', reqacc_cll
435  write(ncheckout_cll,*) 'critacc ', critacc_cll
436  write(ncheckout_cll,*) 'checkacc ', checkacc_cll
437  write(ncheckout_cll,*) 'ErrFlag ', errflag_cll
438  write(ncheckout_cll,*) '---------------------------------------------------------------------------'
439  do i=1,n-1
440  do j=0,3
441  write(ncheckout_cll,fmt3) 'MomVec', j, i, momvec(j,i)
442  end do
443  end do
444  do i=1,binomtable(2,n)
445  write(ncheckout_cll,fmt1) 'MomInv ', i, mominv(i)
446  end do
447  do i=0,n-1
448  write(ncheckout_cll,fmt1) 'masses2', i, masses2(i)
449  end do
450  write(ncheckout_cll,*) '---------------------------------------------------------------------------'
451  write(ncheckout_cll,fmt2) 'COLI:',0,0,0,0,t1(0,0,0,0),0
452  write(ncheckout_cll,fmt2) 'DD :',0,0,0,0,t2(0,0,0,0),0
453  write(ncheckout_cll,fmt2) 'COLI:',n0,n1,n2,n3,t1(n0,n1,n2,n3),r
454  write(ncheckout_cll,fmt2) 'DD :',n0,n1,n2,n3,t2(n0,n1,n2,n3),r
455  if(norm(r).ne.0d0)then
456  write(ncheckout_cll,*) 'diff:', abs(difftn)/norm(r)
457  ratio=abs(difftn)/norm(r)
458  else
459  write(ncheckout_cll,*) 'diff:', 1d50
460  ratio=1d50
461  endif
462  flag=2
463  elseif((flag.eq.2).and.(abs(difftn).gt.ratio*norm(r))) then
464  write(ncheckout_cll,fmt2) 'COLI:',n0,n1,n2,n3,t1(n0,n1,n2,n3),r
465  write(ncheckout_cll,fmt2) 'DD :',n0,n1,n2,n3,t2(n0,n1,n2,n3),r
466  if(norm(r).gt.1d-100)then
467  write(ncheckout_cll,*) 'diff:', abs(difftn)/norm(r)
468  ratio=abs(difftn)/norm(r)
469  else
470  write(ncheckout_cll,*) 'diff:', 1d50
471  ratio=1d50
472  endif
473  elseif ((flag.eq.0).and.(abs(difftn).gt.checkacc_cll*norm(r))) then
474  flag=3
475  end if
476 
477  end do
478  end do
479  end do
480  end do
481 
482  if(flag.eq.2)then
483  write(ncheckout_cll,*) '*************************************************************************'
484  write(ncheckout_cll,*) ' end TNten '
485  write(ncheckout_cll,*)
486  diffcntten_cll(n) = diffcntten_cll(n) + 1
487  if(diffcntten_cll(n).eq.maxcheck_cll(n)) then
488  write(ncheckout_cll,'((A),I4)') ' Further output for differences in TNten functions suppressed for N =', n
489  write(ncheckout_cll,*)
490  endif
491  else if (flag.eq.3) then
492  diffcntten_cll(n) = diffcntten_cll(n) + 1
493  endif
494 
495  end subroutine checktensors_cll
496 
497 
498 
499 
500 
501  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
502  ! subroutine CheckTensorsList_cll(T1,T2,MomVec,MomInv,masses2,N,rmax,Tdiff)
503  !
504  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
505 
506  subroutine checktensorslist_cll(T1,T2,MomVec,MomInv,masses2,norm,N,rmax,Tdiff)
507 
508  integer, intent(in) :: N, rmax
509 ! double complex, intent(in) :: T1(0:RtS(rmax)), T2(0:RtS(rmax))
510  double complex, intent(in) :: T1(RtS(rmax)), T2(RtS(rmax))
511  double complex, intent(in) :: MomVec(0:3,N-1)
512  double complex, intent(in) :: MomInv(BinomTable(2,N)), masses2(0:N-1)
513  double precision, intent(in) :: norm(0:rmax)
514  double precision, intent(out) :: Tdiff(0:rmax)
515  integer :: r,ind,i,j,flag,n0,n1,n2,n3
516  double complex :: diffTN
517  double precision :: ratio
518  character(len=*),parameter :: fmt1 = "(A8,'(',i2,') = dcmplx(',d25.18,' , ',d25.18,' )')"
519  character(len=*),parameter :: fmt2 = &
520  "(A6,' TNten(',i1,',',i1,',',i1,',',i1,') = (',d23.16,' , ',d23.16,' ), r=',i2)"
521  character(len=*),parameter :: fmt3 = &
522  "(A7,'(',i1,',',i2,') = dcmplx(',d25.18,' , ',d25.18,' )')"
523 
524  checkcntten_cll(n) = checkcntten_cll(n) + 1
525 
526 ! data CheckCntTN /0/
527 
528  flag=1
529  if(diffcntten_cll(n).ge.maxcheck_cll(n)) flag=0
530  if(ncheckout_cll.eq.-1) flag=0
531 
532  ratio=0d0
533 
534  tdiff=0d0
535  do r=0,rmax
536  do ind=rts(r-1)+1,rts(r)
537 
538  n0 = lorindtab(0,ind)
539  n1 = lorindtab(1,ind)
540  n2 = lorindtab(2,ind)
541  n3 = lorindtab(3,ind)
542 
543  difftn = t1(ind)-t2(ind)
544  tdiff(r) = max(tdiff(r),abs(difftn)/norm(r))
545 
546  if ((abs(dreal(difftn)).gt.checkacc_cll*norm(r).or.abs(dimag(difftn)).gt.checkacc_cll*norm(r)) .and.(flag.eq.1)) then
547  write(ncheckout_cll,*) '***************************************************************************'
548  write(ncheckout_cll,'(A15,I2,A16,I10)') 'TNten with N =',n,' difference NO.', diffcntten_cll(n)+1
549  write(ncheckout_cll,*) 'COLI and DD do not agree! checkacc =', checkacc_cll
550  write(ncheckout_cll,'(A24,I2,A10,I2,A4,I2)') 'TNten integral with N =', n, ' and rank ', r,' of ',rmax
551  write(ncheckout_cll,*) '---------------------------------------------------------------------------'
552  write(ncheckout_cll,*) 'GLOBAL PARAMETERS:'
553  write(ncheckout_cll,*) 'mode ', mode_cll
554  write(ncheckout_cll,*) 'muUV2 ', muuv2_cll
555  write(ncheckout_cll,*) 'muIR2 ', muir2_cll
556  write(ncheckout_cll,*) 'deltaUV ', deltauv_cll
557  write(ncheckout_cll,*) 'deltaIR1 ', deltair1_cll
558  write(ncheckout_cll,*) 'deltaIR2 ', deltair2_cll
559  write(ncheckout_cll,*) 'nminf ', nminf_cll
560  do i=1,nminf_cll
561  write(ncheckout_cll,*) 'minf2 ', i, minf2_cll(i)
562  end do
563  write(ncheckout_cll,*) 'dprec ', dprec_cll
564  write(ncheckout_cll,*) 'reqacc ', reqacc_cll
565  write(ncheckout_cll,*) 'critacc ', critacc_cll
566  write(ncheckout_cll,*) 'checkacc ', checkacc_cll
567  write(ncheckout_cll,*) 'ErrFlag ', errflag_cll
568  write(ncheckout_cll,*) '---------------------------------------------------------------------------'
569  do i=1,n-1
570  do j=0,3
571  write(ncheckout_cll,fmt3) 'MomVec', j, i, momvec(j,i)
572  end do
573  end do
574  do i=1,binomtable(2,n)
575  write(ncheckout_cll,fmt1) 'MomInv ', i, mominv(i)
576  end do
577  do i=0,n-1
578  write(ncheckout_cll,fmt1) 'masses2', i, masses2(i)
579  end do
580  write(ncheckout_cll,*) '---------------------------------------------------------------------------'
581  write(ncheckout_cll,fmt2) 'COLI:',0,0,0,0,t1(1),0
582  write(ncheckout_cll,fmt2) 'DD :',0,0,0,0,t2(1),0
583  write(ncheckout_cll,fmt2) 'COLI:',n0,n1,n2,n3,t1(ind),r
584  write(ncheckout_cll,fmt2) 'DD :',n0,n1,n2,n3,t2(ind),r
585  if(norm(r).ne.0d0)then
586  write(ncheckout_cll,*) 'diff:', abs(difftn)/norm(r)
587  ratio=abs(difftn)/norm(r)
588  else
589  write(ncheckout_cll,*) 'diff:', 1d50
590  ratio=1d50
591  endif
592  flag=2
593  elseif((flag.eq.2).and.(abs(difftn).gt.ratio*norm(r))) then
594  write(ncheckout_cll,fmt2) 'COLI:',n0,n1,n2,n3,t1(ind),r
595  write(ncheckout_cll,fmt2) 'DD :',n0,n1,n2,n3,t2(ind),r
596  if(norm(r).gt.1d-100)then
597  write(ncheckout_cll,*) 'diff:', abs(difftn)/norm(r)
598  ratio=abs(difftn)/norm(r)
599  else
600  write(ncheckout_cll,*) 'diff:', 1d50
601  ratio=1d50
602  endif
603  elseif ((flag.eq.0).and.(abs(difftn).gt.checkacc_cll*norm(r))) then
604  flag=3
605  end if
606 
607  end do
608  end do
609 
610  if(flag.eq.2)then
611  write(ncheckout_cll,*) '*************************************************************************'
612  write(ncheckout_cll,*) ' end TNten '
613  write(ncheckout_cll,*)
614  diffcntten_cll(n) = diffcntten_cll(n) + 1
615  if(diffcntten_cll(n).eq.maxcheck_cll(n)) then
616  write(ncheckout_cll,*) ' Further output for differences in TNten functions suppressed for N =', n
617  write(ncheckout_cll,*)
618  endif
619  else if (flag.eq.3) then
620  diffcntten_cll(n) = diffcntten_cll(n) + 1
621  endif
622 
623  end subroutine checktensorslist_cll
624 
625 
626 
627 
628 
629  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
630  ! subroutine CheckTenA_cll(T1,T2,mass2,norm,rmax,Tdiff)
631  !
632  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
633 
634  subroutine checktena_cll(T1,T2,masses2,norm,rmax,Tdiff)
635 
636  integer, intent(in) :: rmax
637  double complex, intent(in) :: T1(0:rmax,0:rmax,0:rmax,0:rmax)
638  double complex, intent(in) :: T2(0:rmax,0:rmax,0:rmax,0:rmax)
639  double complex, intent(in) :: masses2(0:0)
640  double precision, intent(in) :: norm(0:rmax)
641  double precision, intent(out) :: Tdiff(0:rmax)
642  integer :: r,n0,n1,n2,n3,i,j,flag
643  double complex :: diffTN
644  double precision :: ratio
645  character(len=*),parameter :: fmt1 = "(A8,'(',i2,') = dcmplx(',d25.18,',',d25.18,' )')"
646  character(len=*),parameter :: fmt2 = &
647  "(A6,' TNten(',i1,',',i1,',',i1,',',i1,') = (',d23.16,' , ',d23.16,' ), r=',i2)"
648  character(len=*),parameter :: fmt3 = &
649  "(A7,'(',i1,',',i2,') = dcmplx(',d25.18,' , ',d25.18,' )')"
650 
651  checkcntten_cll(1) = checkcntten_cll(1) + 1
652 
653 ! data CheckCntTN /0/
654 
655  flag=1
656  if(diffcntten_cll(1).ge.maxcheck_cll(1)) flag=0
657  if(ncheckout_cll.eq.-1) flag=0
658 
659  ratio=0d0
660 
661  tdiff=0d0
662  do r=0,rmax
663  do n0=0,r
664  do n1=0,r-n0
665  do n2=0,r-n0-n1
666  n3=r-n0-n1-n2
667 
668  difftn = t1(n0,n1,n2,n3)-t2(n0,n1,n2,n3)
669  tdiff(r) = max(tdiff(r),abs(difftn)/norm(r))
670  if ((abs(dreal(difftn)).gt.checkacc_cll*norm(r).or.abs(dimag(difftn)).gt.checkacc_cll*norm(r)) .and.(flag.eq.1)) then
671  write(ncheckout_cll,*) '***************************************************************************'
672  write(ncheckout_cll,'(A15,I2,A16,I10)') 'TNten with N =',1,' difference NO.', diffcntten_cll(1)+1
673  write(ncheckout_cll,*) 'COLI and DD do not agree! checkacc =', checkacc_cll
674  write(ncheckout_cll,'(A24,I2,A10,I2,A4,I2)') 'TNten integral with N =', 1, ' and rank ', r,' of ',rmax
675  write(ncheckout_cll,*) '---------------------------------------------------------------------------'
676  write(ncheckout_cll,*) 'GLOBAL PARAMETERS:'
677  write(ncheckout_cll,*) 'mode ', mode_cll
678  write(ncheckout_cll,*) 'muUV2 ', muuv2_cll
679  write(ncheckout_cll,*) 'muIR2 ', muir2_cll
680  write(ncheckout_cll,*) 'deltaUV ', deltauv_cll
681  write(ncheckout_cll,*) 'deltaIR1 ', deltair1_cll
682  write(ncheckout_cll,*) 'deltaIR2 ', deltair2_cll
683  write(ncheckout_cll,*) 'nminf ', nminf_cll
684  do i=1,nminf_cll
685  write(ncheckout_cll,*) 'minf2 ', i, minf2_cll(i)
686  end do
687  write(ncheckout_cll,*) 'dprec ', dprec_cll
688  write(ncheckout_cll,*) 'reqacc ', reqacc_cll
689  write(ncheckout_cll,*) 'critacc ', critacc_cll
690  write(ncheckout_cll,*) 'checkacc ', checkacc_cll
691  write(ncheckout_cll,*) 'ErrFlag ', errflag_cll
692  write(ncheckout_cll,*) '---------------------------------------------------------------------------'
693  write(ncheckout_cll,fmt1) 'masses2', 0, masses2(0)
694  write(ncheckout_cll,*) '---------------------------------------------------------------------------'
695  write(ncheckout_cll,fmt2) 'COLI:',0,0,0,0,t1(0,0,0,0),0
696  write(ncheckout_cll,fmt2) 'DD :',0,0,0,0,t2(0,0,0,0),0
697  write(ncheckout_cll,fmt2) 'COLI:',n0,n1,n2,n3,t1(n0,n1,n2,n3),r
698  write(ncheckout_cll,fmt2) 'DD :',n0,n1,n2,n3,t2(n0,n1,n2,n3),r
699  if(norm(r).ne.0d0)then
700  write(ncheckout_cll,*) 'diff:', abs(difftn)/norm(r)
701  ratio=abs(difftn)/norm(r)
702  else
703  write(ncheckout_cll,*) 'diff:', 1d50
704  ratio=1d50
705  endif
706  flag=2
707  elseif((flag.eq.2).and.(abs(difftn).gt.ratio*norm(r))) then
708  write(ncheckout_cll,fmt2) 'COLI:',n0,n1,n2,n3,t1(n0,n1,n2,n3),r
709  write(ncheckout_cll,fmt2) 'DD :',n0,n1,n2,n3,t2(n0,n1,n2,n3),r
710  if(norm(r).gt.1d-100)then
711  write(ncheckout_cll,*) 'diff:', abs(difftn)/norm(r)
712  ratio=abs(difftn)/norm(r)
713  else
714  write(ncheckout_cll,*) 'diff:', 1d50
715  ratio=1d50
716  endif
717  elseif ((flag.eq.0).and.(abs(difftn).gt.checkacc_cll*norm(r))) then
718  flag=3
719  end if
720 
721  end do
722  end do
723  end do
724  end do
725 
726  if(flag.eq.2)then
727  write(ncheckout_cll,*) '*************************************************************************'
728  write(ncheckout_cll,*) ' end TNten '
729  write(ncheckout_cll,*)
730  diffcntten_cll(1) = diffcntten_cll(1) + 1
731  if(diffcntten_cll(1).eq.maxcheck_cll(1)) then
732  write(ncheckout_cll,*) ' Further output for differences in TNten functions suppressed for N =', 1
733  write(ncheckout_cll,*)
734  endif
735  else if (flag.eq.3) then
736  diffcntten_cll(1) = diffcntten_cll(1) + 1
737  endif
738 
739  end subroutine checktena_cll
740 
741 
742 
743 
744 
745  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
746  ! subroutine CheckTensorsAList_cll(T1,T2,masses2,norm,rmax,Tdiff)
747  !
748  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
749 
750  subroutine checktenalist_cll(T1,T2,masses2,norm,rmax,Tdiff)
751 
752  integer, intent(in) :: rmax
753 ! double complex, intent(in) :: T1(0:RtS(rmax)), T2(0:RtS(rmax))
754  double complex, intent(in) :: T1(RtS(rmax)), T2(RtS(rmax))
755  double complex, intent(in) :: masses2(0:0)
756  double precision, intent(in) :: norm(0:rmax)
757  double precision, intent(out) :: Tdiff(0:rmax)
758  integer :: r,ind,i,j,flag,n0,n1,n2,n3
759  double complex :: diffTN
760  double precision :: ratio
761  character(len=*),parameter :: fmt1 = "(A8,'(',i2,') = dcmplx(',d25.18,',',d25.18,' )')"
762  character(len=*),parameter :: fmt2 = &
763  "(A6,' TNten(',i1,',',i1,',',i1,',',i1,') = (',d23.16,' , ',d23.16,' ), r=',i2)"
764  character(len=*),parameter :: fmt3 = &
765  "(A7,'(',i1,',',i2,') = dcmplx(',d25.18,' , ',d25.18,' )')"
766 
767  checkcntten_cll(1) = checkcntten_cll(1) + 1
768 
769 ! data CheckCntTN /0/
770 
771  flag=1
772  if(diffcntten_cll(1).ge.maxcheck_cll(1)) flag=0
773  if(ncheckout_cll.eq.-1) flag=0
774 
775  ratio=0d0
776 
777  tdiff=0d0
778  do r=0,rmax
779  do ind=rts(r-1)+1,rts(r)
780 
781  n0 = lorindtab(0,ind)
782  n1 = lorindtab(1,ind)
783  n2 = lorindtab(2,ind)
784  n3 = lorindtab(3,ind)
785  difftn = t1(ind)-t2(ind)
786  tdiff(r) = max(tdiff(r),abs(difftn)/norm(r))
787 
788  if ((abs(dreal(difftn)).gt.checkacc_cll*norm(r).or.abs(dimag(difftn)).gt.checkacc_cll*norm(r)) .and.(flag.eq.1)) then
789  write(ncheckout_cll,*) '***************************************************************************'
790  write(ncheckout_cll,'(A15,I2,A16,I10)') 'TNten with N =',1,' difference NO.', diffcntten_cll(1)+1
791  write(ncheckout_cll,*) 'COLI and DD do not agree! checkacc =', checkacc_cll
792  write(ncheckout_cll,'(A24,I2,A10,I2,A4,I2)') 'TNten integral with N =', 1, ' and rank ', r,' of ',rmax
793  write(ncheckout_cll,*) '---------------------------------------------------------------------------'
794  write(ncheckout_cll,*) 'GLOBAL PARAMETERS:'
795  write(ncheckout_cll,*) 'mode ', mode_cll
796  write(ncheckout_cll,*) 'muUV2 ', muuv2_cll
797  write(ncheckout_cll,*) 'muIR2 ', muir2_cll
798  write(ncheckout_cll,*) 'deltaUV ', deltauv_cll
799  write(ncheckout_cll,*) 'deltaIR1 ', deltair1_cll
800  write(ncheckout_cll,*) 'deltaIR2 ', deltair2_cll
801  write(ncheckout_cll,*) 'nminf ', nminf_cll
802  do i=1,nminf_cll
803  write(ncheckout_cll,*) 'minf2 ', i, minf2_cll(i)
804  end do
805  write(ncheckout_cll,*) 'dprec ', dprec_cll
806  write(ncheckout_cll,*) 'reqacc ', reqacc_cll
807  write(ncheckout_cll,*) 'critacc ', critacc_cll
808  write(ncheckout_cll,*) 'checkacc ', checkacc_cll
809  write(ncheckout_cll,*) 'ErrFlag ', errflag_cll
810  write(ncheckout_cll,*) '---------------------------------------------------------------------------'
811  write(ncheckout_cll,fmt1) 'masses2', 0, masses2(0)
812  write(ncheckout_cll,*) '---------------------------------------------------------------------------'
813  write(ncheckout_cll,fmt2) 'COLI:',0,0,0,0,t1(1),0
814  write(ncheckout_cll,fmt2) 'DD :',0,0,0,0,t2(1),0
815  write(ncheckout_cll,fmt2) 'COLI:',n0,n1,n2,n3,t1(ind),r
816  write(ncheckout_cll,fmt2) 'DD :',n0,n1,n2,n3,t2(ind),r
817  if(norm(r).ne.0d0)then
818  write(ncheckout_cll,*) 'diff:', abs(difftn)/norm(r)
819  ratio=abs(difftn)/norm(r)
820  else
821  write(ncheckout_cll,*) 'diff:', 1d50
822  ratio=1d50
823  endif
824  flag=2
825  elseif((flag.eq.2).and.(abs(difftn).gt.ratio*norm(r))) then
826  write(ncheckout_cll,fmt2) 'COLI:',n0,n1,n2,n3,t1(ind),r
827  write(ncheckout_cll,fmt2) 'DD :',n0,n1,n2,n3,t2(ind),r
828  if(norm(r).gt.1d-100)then
829  write(ncheckout_cll,*) 'diff:', abs(difftn)/norm(r)
830  ratio=abs(difftn)/norm(r)
831  else
832  write(ncheckout_cll,*) 'diff:', 1d50
833  ratio=1d50
834  endif
835  elseif ((flag.eq.0).and.(abs(difftn).gt.checkacc_cll*norm(r))) then
836  flag=3
837  end if
838 
839  end do
840  end do
841 
842  if(flag.eq.2)then
843  write(ncheckout_cll,*) '*************************************************************************'
844  write(ncheckout_cll,*) ' end TNten '
845  write(ncheckout_cll,*)
846  diffcntten_cll(1) = diffcntten_cll(1) + 1
847  if(diffcntten_cll(1).eq.maxcheck_cll(1)) then
848  write(ncheckout_cll,*) ' Further output for differences in TNten functions suppressed for N =', 1
849  write(ncheckout_cll,*)
850  endif
851  else if (flag.eq.3) then
852  diffcntten_cll(1) = diffcntten_cll(1) + 1
853  endif
854 
855  end subroutine checktenalist_cll
856 
857 
858 end module inittensors
inittensors::time3
real time3
Definition: InitTensors.F90:40
inittensors::indspiprod
integer, dimension(:,:,:,:), allocatable indspiprod
Definition: InitTensors.F90:38
inittensors::checktensorslist_cll
subroutine checktensorslist_cll(T1, T2, MomVec, MomInv, masses2, norm, N, rmax, Tdiff)
Definition: InitTensors.F90:507
inittensors::setlorindtab
subroutine setlorindtab(rmax)
Definition: InitTensors.F90:98
inittensors::rts
integer, dimension(:), allocatable rts
Definition: InitTensors.F90:39
inittensors::stin
real stin
Definition: InitTensors.F90:40
inittensors
Definition: InitTensors.F90:25
inittensors::calcindspiprod
integer function, dimension(:,:,:), allocatable calcindspiprod(Nm1, r)
Definition: InitTensors.F90:302
inittensors::cftab
integer, dimension(:,:), allocatable cftab
Definition: InitTensors.F90:37
inittensors::time1
real time1
Definition: InitTensors.F90:40
inittensors::switchoncalcuv_ten
subroutine switchoncalcuv_ten()
Definition: InitTensors.F90:364
inittensors::checktensors_cll
subroutine checktensors_cll(T1, T2, MomVec, MomInv, masses2, norm, N, rmax, Tdiff)
Definition: InitTensors.F90:379
inittensors::setindspiprod
subroutine setindspiprod(Nm1max, rmax)
Definition: InitTensors.F90:261
inittensors::switchoffcalcuv_ten
subroutine switchoffcalcuv_ten()
Definition: InitTensors.F90:349
inittensors::calcuv_ten
logical calcuv_ten
Definition: InitTensors.F90:41
inittensors::setaddgtab
subroutine setaddgtab(rmax)
Definition: InitTensors.F90:171
inittensors::setrts
subroutine setrts(rmax)
Definition: InitTensors.F90:74
inittensors::checktena_cll
subroutine checktena_cll(T1, T2, masses2, norm, rmax, Tdiff)
Definition: InitTensors.F90:635
inittensors::addgtab
integer, dimension(:,:), allocatable addgtab
Definition: InitTensors.F90:37
inittensors::stout
real stout
Definition: InitTensors.F90:40
inittensors::setcftab
subroutine setcftab(rmax)
Definition: InitTensors.F90:220
inittensors::lorindtab
integer, dimension(:,:), allocatable lorindtab
Definition: InitTensors.F90:37
inittensors::checktenalist_cll
subroutine checktenalist_cll(T1, T2, masses2, norm, rmax, Tdiff)
Definition: InitTensors.F90:751
inittensors::setaddindtab
subroutine setaddindtab(rmax)
Definition: InitTensors.F90:136
inittensors::addindtab
integer, dimension(:,:), allocatable addindtab
Definition: InitTensors.F90:37
inittensors::time2
real time2
Definition: InitTensors.F90:40
inittensors::init_tables2
subroutine init_tables2(Nm1, rmax)
Definition: InitTensors.F90:52