JHUGen MELA  JHUGen v7.5.6, MELA v2.4.2
Matrix element calculations as used in JHUGen.
demo.f90
Go to the documentation of this file.
1 !!
2 !! File demo.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 program demo
11 
12  use collier
13 
14  implicit none
15 
16  double complex :: m02,m12,m22,m32,m42,m52
17  double complex :: p10,p21,p32,p43,p54,p50,p20,p31
18  double complex :: p42,p53,p40,p51,p30,p41,p52
19  double complex :: p1vec(0:3),p2vec(0:3),p3vec(0:3),p4vec(0:3),p5vec(0:3)
20  double complex :: mominv(15), masses2(0:5), momvec(0:3,5)
21  double complex, allocatable :: acoeff(:),acoeffuv(:)
22  double complex, allocatable :: bcoeff(:,:),bcoeffuv(:,:)
23  double complex, allocatable :: ccoeff(:,:,:),ccoeffuv(:,:,:)
24  double complex, allocatable :: dcoeff(:,:,:,:),dcoeffuv(:,:,:,:)
25  double complex, allocatable :: ecoeff(:,:,:,:,:),ecoeffuv(:,:,:,:,:)
26  double complex, allocatable :: fcoeff(:,:,:,:,:,:),fcoeffuv(:,:,:,:,:,:)
27  double complex, allocatable :: cten(:,:,:,:),ctenuv(:,:,:,:)
28  double complex, allocatable :: dten(:,:,:,:),dtenuv(:,:,:,:)
29  double complex, allocatable :: eten(:,:,:,:),etenuv(:,:,:,:)
30  double complex, allocatable :: ften(:,:,:,:),ftenuv(:,:,:,:)
31  double complex, allocatable :: dbcoeff(:,:),dbcoeffuv(:,:)
32  double complex, allocatable :: tncoeff(:),tncoeffuv(:)
33  double complex, allocatable :: tnten(:),tntenuv(:)
34  double complex :: a0,b0,c0,d0,e0,f0,a0uv,b0uv,db0,db1,db00,db00uv
35  double complex, parameter :: c0uv=(0d0,0d0),d0uv=(0d0,0d0)
36  double complex, parameter :: e0uv=(0d0,0d0),f0uv=(0d0,0d0)
37  double precision, allocatable :: aerr(:),berr(:),cerr(:),derr(:)
38  double precision, allocatable :: eerr(:),ferr(:),dberr(:)
39  double precision :: muuv2,muir2,deltauv,deltair1,deltair2
40  integer :: rank,nmax,n,rmax
41  integer :: mode,casen,casei
42 
43 ! Nmax = maximal degree N of N-point function
44  nmax = 6
45  rmax = 6
46 
47 ! Initialize COLLIER
48 ! the default directory for diagnostic messages is "output_cll"
49  call init_cll(nmax,rmax,"output_cll")
50 
51 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
52 ! for checks
53 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
54  call initcachesystem_cll(1,nmax)
55  call initevent_cll
56 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
57 
58 ! uncomment to redefine regularization parameters of COLLIER
59 ! (after initialization)
60 ! MuUV2 = dcmplx(1d2)
61 ! call SetMuIR2_cll(MuUV2)
62 ! MuIR2 = dcmplx(1d2)
63 ! call SetMuIR2_cll(MuIR2)
64 ! DeltaUV = dcmplx(1d0)
65 ! call SetdeltaUV_cll(DeltaUV)
66 ! DeltaIR1 = dcmplx(1d0)
67 ! DeltaIR2 = dcmplx(1d0)
68 ! call SetdeltaIR_cll(DeltaIR1,DeltaIR2)
69 
70 ! uncomment to choose DD instead of COLI library
71 ! call SetMode_cll(2)
72 ! or to perform comparison between COLI and DD
73 ! call setMode_cll(3)
74 
75 ! the specific examples start here
76 
77 ! choose COLLIER mode
78  write(*,*)
79  write(*,*) 'Choose mode od COLLIER :'
80  write (*,*) ' 1) use COLI branch'
81  write (*,*) ' 2) use DD branch'
82  write (*,*) ' 3) use both branches and compare'
83  read(*,*) mode
84 
85 
86  select case(mode)
87 
88  case(1)
89 
90  case(2)
91 
92  case(3)
93 
94  case default
95 
96  write(*,'(/(a))') ' Input does not correspond to available mode'
97 
98  end select
99 
100  call setmode_cll(mode)
101 
102  write(*,*)
103  write(*,*) 'Choose type of function:'
104  write (*,*) ' 1) 1-point function'
105  write (*,*) ' 2) 2-point function'
106  write (*,*) ' 3) 3-point function'
107  write (*,*) ' 4) 4-point function'
108  write (*,*) ' 5) 5-point function'
109  write (*,*) ' 6) 6-point function'
110  write (*,*) '12) derivative of 2-point function'
111  read(*,*) casen
112 
113 
114  select case(casen)
115 
116 ! 1-point function calls
117  case(1)
118 
119  write(*,*)
120  write(*,*) 'Choose test case for 1-point function:'
121  write (*,*) ' 1) regular scalar 1-point function'
122  write (*,*) '10) 1-point coefficients up to rank 2'
123  read(*,*) casei
124 
125  select case(casei)
126 
127 ! example 1 1
128 ! call of regular scalar 1-point function
129  case(1)
130 
131  write(*,'(/(a)/)') ' Calculating regular scalar 1-point function'
132 
133  m02 = 4d2
134 
135  masses2(0) = m02
136 
137 ! A0 contains the results for the scalar 2-point function
138  call a0_cll(a0,m02)
139 
140 ! writes the result to 'demo_1point_example01.dat'
141  call writeresulta0(casen,casei,a0,masses2(0:0))
142 
143 
144 
145 ! example 1 10
146 ! call of tensor 1-point function: coefficients up to rank = 3
147  case (10)
148 
149  write(*,'(/(a))') ' Calculating 1-point coefficients'
150 
151  m02 = 4d2
152 
153  masses2(0:0) = m02
154 
155 ! allocate fields for tensor coefficients
156  rank = 2
157  allocate(acoeff(0:rank/2))
158  allocate(acoeffuv(0:rank/2))
159  allocate(aerr(0:rank))
160 
161 
162 ! Acoeff(n0) contains the results for the corresponding tensor coefficent
163 ! Acoeffuv(n0) contains the results for the UV-singular part of the tensor coefficent
164 ! (coefficient of 2/(4-D)=1/eps term)
165 ! Aerr(r) is optional and contains an absolute estimate of the error of the
166 ! tensor coefficients of rank r (not involving metric tensors)
167  call a_cll(acoeff,acoeffuv,m02,rank,aerr)
168 
169 ! equivalent call
170  n = 1
171  allocate(tncoeff(getnc_cll(n,rank)))
172  allocate(tncoeffuv(getnc_cll(n,rank)))
173  call tn_cll(tncoeff,tncoeffuv,masses2(0:2),n,rank,aerr)
174 
175 ! writes the result to 'demo_1point_example10.dat'
176  call writeresulta(casen,casei,acoeff,acoeffuv,masses2(0:0),rank,aerr)
177 
178  deallocate(acoeff,acoeffuv,aerr,tncoeff,tncoeffuv)
179 
180  case default
181 
182  write(*,'(/(a))') ' Input does not correspond to predefined sample'
183 
184  end select
185 
186 
187 ! 2-point function calls
188  case(2)
189 
190  write(*,*)
191  write(*,*) 'Choose test case for 2-point function:'
192  write (*,*) ' 1) regular scalar 2-point function'
193  write (*,*) ' 2) scalar 2-point function with small masses'
194  write (*,*) '10) tensor 2-point function: coefficients up to rank 3'
195  read(*,*) casei
196 
197  select case(casei)
198 
199 ! example 2 1
200 ! call of regular scalar 2-point function
201  case(1)
202 
203  write(*,'(/(a)/)') ' Calculating regular scalar 2-point function'
204 
205  p10 = 1d3
206  m02 = 4d2
207  m12 = 1d2
208 
209  mominv(1:1) = p10
210  masses2(0:1) = (/m02,m12/)
211 
212 ! B0 contains the results for the scalar 2-point function
213  call b0_cll(b0,p10,m02,m12)
214 
215 ! equivalent call
216  call b0_cll(b0,mominv(1:1),masses2(0:1))
217 
218 ! writes the result to 'demo_2point_example01.dat'
219  call writeresultb0(casen,casei,b0,mominv(1:1),masses2(0:1))
220 
221 
222 ! example 2 2
223 ! call of regular scalar 2-point function with small masses
224  case(2)
225 
226  write(*,'(/(a)/)') &
227  ' Calculating scalar 2-point function with small masses'
228 
229  p10 = 4d2
230  m02 = 1d-4
231  m12 = 1d-4
232 
233 ! clears list of regulator masses
234  call clearminf2_cll
235 ! adds m12 to list of regulator masses
236 ! these are treated as infinitesimal apart from their appearance
237 ! in IR-singular logarithms, where the actual values are taken
238  call addminf2_cll(m12)
239 
240  mominv(1:1) = p10
241  masses2(0:1) = (/m02,m12/)
242 
243 ! B0 contains the results for the scalar 2-point function
244  call b0_cll(b0,p10,m02,m12)
245 
246 ! equivalent call
247  call b0_cll(b0,mominv(1:1),masses2(0:1))
248 
249 ! writes the result to 'demo_2point_example02.dat'
250  call writeresultb0(casen,casei,b0,mominv(1:1),masses2(0:1))
251 
252 
253 ! example 2 10
254 ! call of tensor 2-point function: coefficients up to rank = 3
255  case (10)
256 
257  write(*,'(/(a))') ' Calculating 2-point coefficients'
258 
259  p10 = 1d3
260  m02 = 4d2
261  m12 = 1d2
262 
263  mominv(1:1) = (/p10/)
264  masses2(0:1) = (/m02,m12/)
265 
266 ! allocate fields for tensor coefficients
267  rank = 3
268  allocate(bcoeff(0:rank/2,0:rank))
269  allocate(bcoeffuv(0:rank/2,0:rank))
270  allocate(berr(0:rank))
271 
272 
273 ! Bcoeff(n0,n1) contains the results for the corresponding tensor coefficent
274 ! Bcoeffuv(n0,n1) contains the results for the UV-singular part of the tensor coefficent
275 ! (coefficient of 2/(4-D)=1/eps term)
276 ! Berr(r) is optional and contains an absolute estimate of the error of the
277 ! tensor coefficients of rank r (not involving metric tensors)
278  call b_cll(bcoeff,bcoeffuv,p10,m02,m12,rank,berr)
279 
280 ! equivalent call
281  call b_cll(bcoeff,bcoeffuv,mominv(1:1),masses2(0:1),rank,berr)
282 
283 ! equivalent calls, but all coefficients are in 1-dimensional array TNcoeff
284  n = 2
285  allocate(tncoeff(getnc_cll(n,rank)))
286  allocate(tncoeffuv(getnc_cll(n,rank)))
287  call b_cll(tncoeff,tncoeffuv,p10,m02,m12,rank,berr)
288 ! or
289  call b_cll(tncoeff,tncoeffuv,mominv(1:1),masses2(0:1),rank,berr)
290 ! or
291  call tn_cll(tncoeff,tncoeffuv,mominv(1:1),masses2(0:1),n,rank,berr)
292 
293 ! writes the result to 'demo_2point_example10.dat'
294  call writeresultb(casen,casei,bcoeff,bcoeffuv,mominv(1:1),masses2(0:1),rank,berr)
295 
296  deallocate(bcoeff,bcoeffuv,berr,tncoeff,tncoeffuv)
297 
298  case default
299 
300  write(*,'(/(a))') ' Input does not correspond to predefined sample'
301 
302  end select
303 
304 ! 3-point function calls
305  case(3)
306 
307  write(*,*)
308  write(*,*) 'Choose test case for 3-point function:'
309  write (*,*) ' 1) regular scalar 3-point function'
310  write (*,*) ' 2) IR-singular scalar 3-point function'
311  write (*,*) ' 3) Mass-singular scalar 3-point function in dimensional regularization'
312  write (*,*) ' 4) Mass-singular scalar 3-point function in mass regularization'
313  write (*,*) '10) tensor 3-point function: coefficients up to rank 3'
314  write (*,*) '20) tensor 3-point function: tensor components up to rank 3'
315  read(*,*) casei
316 
317  select case(casei)
318 
319 ! example 3 1
320 ! call of regular scalar 3-point function
321  case(1)
322 
323  write(*,'(/(a)/)') ' Calculating regular scalar 3-point function'
324 
325  p10 = 2d2
326  p21 = 1d3
327  p20 = 2d2
328  m02 = 1d2
329  m12 = 1d2
330  m22 = 1d2
331 
332  mominv(1:3) = (/p10,p21,p20/)
333  masses2(0:2) = (/m02,m12,m22/)
334 
335 ! C0 contains the results for the scalar 3-point function
336  call c0_cll(c0,p10,p21,p20,m02,m12,m22)
337 
338 ! equivalent call
339  call c0_cll(c0,mominv(1:3),masses2(0:2))
340 
341 ! writes the result to 'demo_3point_example01.dat'
342  call writeresultc0(casen,casei,c0,mominv(1:3),masses2(0:2))
343 
344 
345 ! example 3 2
346 ! call of IR-singular scalar 3-point function
347  case(2)
348 
349  write(*,'(/(a)/)') ' Calculating IR-singular scalar 3-point function'
350  write(*,'((a))') ' IR singularity described by log(muir2) and deltair'
351 
352  p10 = 1d2
353  p21 = 5d2
354  p20 = 1d2
355  m02 = 0d0
356  m12 = 1d2
357  m22 = 1d2
358 
359  mominv(1:3) = (/p10,p21,p20/)
360  masses2(0:2) = (/m02,m12,m22/)
361 
362 ! C0 contains the results for the scalar 3-point function
363  call c0_cll(c0,p10,p21,p20,m02,m12,m22)
364 
365 ! equivalent call
366  call c0_cll(c0,mominv(1:3),masses2(0:2))
367 
368 ! writes the result to 'demo_3point_example02.dat'
369  call writeresultc0(casen,casei,c0,mominv(1:3),masses2(0:2))
370 
371 
372 ! example 3 3
373 ! call of mass-singular scalar 3-point function in dimensional regularization
374  case(3)
375 
376  write(*,'(/(a))') ' Calculating mass-singular scalar 3-point function'
377  write(*,'((a)/)') ' in dimensional regularization'
378  write(*,'((a))') ' mass-singularity described by log(muir2) and deltair'
379 
380  p10 = 0d0
381  p21 = 5d2
382  p20 = 4d2
383  m02 = 0d0
384  m12 = 0d0
385  m22 = 0d0
386 
387  mominv(1:3) = (/p10,p21,p20/)
388  masses2(0:2) = (/m02,m12,m22/)
389 
390 ! C0 contains the results for the scalar 3-point function
391  call c0_cll(c0,p10,p21,p20,m02,m12,m22)
392 
393 ! equivalent call
394  call c0_cll(c0,mominv(1:3),masses2(0:2))
395 
396 ! writes the result to 'demo_3point_example03.dat'
397  call writeresultc0(casen,casei,c0,mominv(1:3),masses2(0:2))
398 
399 
400 ! example 3 4
401 ! call of mass-singular scalar 3-point function in mass regularization
402  case(4)
403 
404  write(*,'(/(a))') ' Calculating mass-singular scalar 3-point function'
405  write(*,'((a)/)') ' in mass regularization'
406  write(*,'((a)/)') ' mass-singularity described by log(m12)'
407 
408  p10 = 1d-4
409  p21 = 5d2
410  p20 = 4d2
411  m02 = 0d0
412  m12 = 1d-4
413  m22 = 0d0
414 
415 ! clears list of regulator masses
416  call clearminf2_cll
417 ! adds m12 to list of regulator masses
418 ! these are treated as infinitesimal apart from their appearance
419 ! in IR-singular logarithms, where the actual values are taken
420  call addminf2_cll(m12)
421 
422  mominv(1:3) = (/p10,p21,p20/)
423  masses2(0:2) = (/m02,m12,m22/)
424 
425 ! C0 contains the results for the scalar 3-point function
426  call c0_cll(c0,p10,p21,p20,m02,m12,m22)
427 
428 ! equivalent call
429  call c0_cll(c0,mominv(1:3),masses2(0:2))
430 
431 ! writes the result to 'demo_3point_example04.dat'
432  call writeresultc0(casen,casei,c0,mominv(1:3),masses2(0:2))
433 
434 
435 
436 ! example 3 10
437 ! call of tensor 3-point function: coefficients up to rank = 3
438  case (10)
439 
440  write(*,'(/(a))') ' Calculating 3-point coefficients'
441 
442  p10 = 2d2
443  p21 = 1d3
444  p20 = 2d2
445  m02 = 1d2
446  m12 = 1d2
447  m22 = 1d2
448 
449  mominv(1:3) = (/p10,p21,p20/)
450  masses2(0:2) = (/m02,m12,m22/)
451 
452 ! allocate fields for tensor coefficients
453  rank = 3
454  allocate(ccoeff(0:rank/2,0:rank,0:rank))
455  allocate(ccoeffuv(0:rank/2,0:rank,0:rank))
456  allocate(cerr(0:rank))
457 
458 
459 ! Ccoeff(n0,n1,n2) contains the results for the corresponding tensor coefficent
460 ! Ccoeffuv(n0,n1,n2) contains the results for the UV-singular part of the tensor coefficent
461 ! (coefficient of 2/(4-D)=1/eps term)
462 ! Cerr(r) is optional and contains an absolute estimate of the error of the
463 ! tensor coefficients of rank r (not involving metric tensors)
464  call c_cll(ccoeff,ccoeffuv,p10,p21,p20,m02,m12,m22,rank,cerr)
465 
466 ! equivalent call
467  call c_cll(ccoeff,ccoeffuv,mominv(1:3),masses2(0:2),rank,cerr)
468 
469 ! equivalent calls, but all coefficients are in 1-dimensional array TNcoeff
470  n = 3
471  allocate(tncoeff(getnc_cll(n,rank)))
472  allocate(tncoeffuv(getnc_cll(n,rank)))
473  call c_cll(tncoeff,tncoeffuv,p10,p21,p20,m02,m12,m22,rank,cerr)
474 ! or
475  call c_cll(tncoeff,tncoeffuv,mominv(1:3),masses2(0:2),rank,cerr)
476 ! or
477  call tn_cll(tncoeff,tncoeffuv,mominv(1:3),masses2(0:2),n,rank,cerr)
478 
479 ! writes the result to 'demo_3point_example10.dat'
480  call writeresultc(casen,casei,ccoeff,ccoeffuv,mominv(1:3),masses2(0:2),rank,cerr)
481 
482  deallocate(ccoeff,ccoeffuv,cerr,tncoeff,tncoeffuv)
483 
484 ! example 3 20
485 ! call of tensor 3-point function: tensor components up to rank = 3
486  case (20)
487 
488  write(*,'(/(a))') ' Calculating 3-point tensor components'
489 
490  p1vec = (/ 15d0,0d0,0d0,5d0 /)
491  p2vec = (/ -15d0,0d0,0d0,5d0 /)
492  p10 = 2d2
493  p21 = 9d2
494  p20 = 2d2
495  m02 = 1d2
496  m12 = 1d2
497  m22 = 1d2
498 
499  momvec(0:3,1) = p1vec
500  momvec(0:3,2) = p2vec
501  mominv(1:3) = (/p10,p21,p20/)
502  masses2(0:2) = (/m02,m12,m22/)
503 
504 ! allocate fields for tensor coefficients
505  n = 3
506  rank = 3
507  allocate(cten(0:rank,0:rank,0:rank,0:rank))
508  allocate(ctenuv(0:rank,0:rank,0:rank,0:rank))
509  allocate(cerr(0:rank))
510 
511 ! Cten(n0,n1,n2,n3) contains the results for the corresponding tensor components
512 ! Ctenuv(n0,n1,n2,n3) contains the results for the UV-singular part of the tensor component
513 ! (coefficient of 2/(4-D)=1/eps term)
514 ! Cerr(r) is optional and contains an absolute estimate of the error of the
515 ! tensor component of rank r (not involving metric tensors)
516  call cten_cll(cten,ctenuv,p1vec,p2vec,p10,p21,p20,m02,m12,m22,rank,cerr)
517 
518 ! equivalent calls
519  call cten_cll(cten,ctenuv,momvec(0:3,1:2),mominv(1:3),masses2(0:2),rank,cerr)
520 ! or
521  call tnten_cll(cten,ctenuv,momvec(0:3,1:2),mominv(1:3),masses2(0:2),n,rank,cerr)
522 
523 ! equivalent calls, but with all components are in 1-dimensional array TNten
524  allocate(tnten(getnt_cll(rank)))
525  allocate(tntenuv(getnt_cll(rank)))
526  call cten_cll(tnten,tntenuv,p1vec,p2vec,p10,p21,p20,m02,m12,m22,rank,cerr)
527 ! or
528  call cten_cll(tnten,tntenuv,momvec(0:3,1:2),mominv(1:3),masses2(0:2),rank,cerr)
529 ! or
530  call tnten_cll(tnten,tntenuv,momvec(0:3,1:2),mominv(1:3),masses2(0:2),n,rank,cerr)
531 
532 ! writes the result to 'demo_3point_example20.dat'
533  call writeresultcten(casen,casei,cten,ctenuv,momvec(0:3,1:2),mominv(1:3),masses2(0:2),rank,cerr)
534 
535  deallocate(cten,ctenuv,cerr,tnten,tntenuv)
536 
537  case default
538 
539  write(*,'(/(a))') ' Input does not correspond to predefined sample'
540 
541  end select
542 
543 ! scalar 4-point function calls
544  case(4)
545 
546  write(*,*)
547  write(*,*) 'Choose test case for 4-point function:'
548  write (*,*) ' 1) regular scalar 4-point function'
549 ! write (*,*) ' 2) IR-singular scalar 4-point function'
550 ! write (*,*) ' 3) Mass-singular scalar 4-point function in dimensional regularization'
551 ! write (*,*) ' 4) Mass-singular scalar 4-point function in mass regularization'
552  write (*,*) '10) tensor 4-point function: coefficients up to rank 3'
553  write (*,*) '20) tensor 4-point function: tensor components up to rank 3'
554  read(*,*) casei
555 
556  select case(casei)
557 
558 ! example 4 1
559 ! call of regular scalar 4-point function
560  case(1)
561 
562  write(*,'(/(a)/)') ' Calculating regular scalar 4-point function'
563 
564  p10 = 1d3
565  p21 = 2d3
566  p32 = 3d3
567  p30 = 1d2
568  p20 = 1d4
569  p31 = -3d3
570  m02 = 1d2
571  m12 = 1d2
572  m22 = 3d2
573  m32 = 3d2
574 
575  mominv(1:6) = (/p10,p21,p32,p30,p20,p31/)
576  masses2(0:3) = (/m02,m12,m22,m32/)
577 
578 ! D0 contains the results for the scalar 4-point function
579  call d0_cll(d0,p10,p21,p32,p30,p20,p31,m02,m12,m22,m32)
580 
581 ! equivalent call
582  call d0_cll(d0,mominv(1:6),masses2(0:3))
583 
584 ! writes the result to 'demo_4point_example01.dat'
585  call writeresultd0(casen,casei,d0,mominv(1:6),masses2(0:3))
586 
587 
588 ! example 4 10
589 ! call of 4-point coefficients up to rank = 3
590  case (10)
591 
592  write(*,'(/(a))') ' Calculating tensor 4-point coefficients'
593 
594  p10 = 1d3
595  p21 = 2d3
596  p32 = 3d3
597  p30 = 1d2
598  p20 = 1d4
599  p31 = -3d3
600  m02 = 1d2
601  m12 = 1d2
602  m22 = 3d2
603  m32 = 3d2
604 
605  mominv(1:6) = (/p10,p21,p32,p30,p20,p31/)
606  masses2(0:3) = (/m02,m12,m22,m32/)
607 
608 ! allocate fields for tensor coefficients
609  rank = 3
610  allocate(dcoeff(0:rank/2,0:rank,0:rank,0:rank))
611  allocate(dcoeffuv(0:rank/2,0:rank,0:rank,0:rank))
612  allocate(derr(0:rank))
613 
614 
615 ! Dcoeff(n0,n1,n2,n3) contains the results for the corresponding tensor coefficent
616 ! Dcoeffuv(n0,n1,n2,n3) contains the results for the UV-singular part of the tensor coefficent
617 ! (coefficient of 2/(4-D)=1/eps term)
618 ! Derr(r) is optional and contains an absolute estimate of the error of the
619 ! tensor coefficients of rank r (not involving metric tensors)
620  call d_cll(dcoeff,dcoeffuv,p10,p21,p32,p30,p20,p31,m02,m12,m22,m32,rank,derr)
621 
622 ! equivalent call
623  call d_cll(dcoeff,dcoeffuv,mominv(1:6),masses2(0:3),rank,derr)
624 
625 ! equivalent calls, but all coefficients are in 1-dimensional array TNcoeff
626  n = 4
627  allocate(tncoeff(getnc_cll(n,rank)))
628  allocate(tncoeffuv(getnc_cll(n,rank)))
629  call d_cll(tncoeff,tncoeffuv,p10,p21,p32,p30,p20,p31,m02,m12,m22,m32,rank,derr)
630 ! or
631  call d_cll(tncoeff,tncoeffuv,mominv(1:6),masses2(0:3),rank,derr)
632 ! or
633  call tn_cll(tncoeff,tncoeffuv,mominv(1:6),masses2(0:3),n,rank,derr)
634 
635 
636 ! writes the result to 'demo_4point_example10.dat'
637  call writeresultd(casen,casei,dcoeff,dcoeffuv,mominv(1:6),masses2(0:3),rank,derr)
638 
639  deallocate(dcoeff,dcoeffuv,derr,tncoeff,tncoeffuv)
640 
641 ! example 4 20
642 ! call of 4-point tensor components up to rank = 3
643  case (20)
644 
645  write(*,'(/(a))') ' Calculating 4-point tensor components'
646 
647  p1vec = (/ 1d2,0d0,0d0,1d2 /)
648  p2vec = (/ 2d2,0d0,0d0,0d0 /)
649  p3vec = (/ 1d2,5d1,0d0,5d1 /)
650  p10 = 0d0
651  p21 = 0d0
652  p32 = 5d3
653  p30 = 5d3
654  p20 = 4d4
655  p31 = -5d3
656  m02 = 1d2
657  m12 = 3d2
658  m22 = 1d2
659  m32 = 3d2
660 
661  momvec(0:3,1) = p1vec
662  momvec(0:3,2) = p2vec
663  momvec(0:3,3) = p3vec
664  mominv(1:6) = (/p10,p21,p32,p30,p20,p31/)
665  masses2(0:3) = (/m02,m12,m22,m32/)
666 
667 ! allocate fields for tensor coefficients
668  n = 4
669  rank = 3
670  allocate(dten(0:rank,0:rank,0:rank,0:rank))
671  allocate(dtenuv(0:rank,0:rank,0:rank,0:rank))
672  allocate(derr(0:rank))
673 
674 
675 ! Dten(n0,n1,n2,n3) contains the results for the corresponding tensor component
676 ! Dtenuv(n0,n1,n2,n3) contains the results for the UV-singular part of the tensor component
677 ! (coefficient of 2/(4-D)=1/eps term)
678 ! Derr(r) is optional and contains an absolute estimate of the error of the
679 ! tensor components of rank r (not involving metric tensors)
680  call dten_cll(dten,dtenuv,p1vec,p2vec,p3vec,p10,p21,p32,p30,p20,p31,m02,m12,m22,m32,rank,derr)
681 
682 ! equivalent calls
683  call dten_cll(dten,dtenuv,momvec(0:3,1:3),mominv(1:6),masses2(0:3),rank,derr)
684 ! or
685  call tnten_cll(dten,dtenuv,momvec(0:3,1:3),mominv(1:6),masses2(0:3),n,rank,derr)
686 
687 ! equivalent calls, but with all components are in 1-dimensional array TNten
688  allocate(tnten(getnt_cll(rank)))
689  allocate(tntenuv(getnt_cll(rank)))
690  call dten_cll(tnten,tntenuv,p1vec,p2vec,p3vec,p10,p21,p32,p30,p20,p31,m02,m12,m22,m32,rank,derr)
691 ! or
692  call dten_cll(tnten,tntenuv,momvec(0:3,1:3),mominv(1:6),masses2(0:3),rank,derr)
693 ! or
694  call tnten_cll(tnten,tntenuv,momvec(0:3,1:3),mominv(1:6),masses2(0:3),n,rank,derr)
695 
696 ! writes the result to 'demo_4point_example20.dat'
697  call writeresultdten(casen,casei,dten,dtenuv,momvec(0:3,1:3),mominv(1:6),masses2(0:3),rank,derr)
698 
699  deallocate(dten,dtenuv,derr,tnten,tntenuv)
700 
701  case default
702 
703  write(*,'(/(a))') ' Input does not correspond to predefined sample'
704 
705  end select
706 
707 ! 5-point function calls
708  case(5)
709 
710  write(*,*)
711  write(*,*) 'Choose test case for 5-point function:'
712  write (*,*) '10) tensor 5-point function: coefficients up to rank 3'
713  write (*,*) '20) tensor 5-point function: tensor components of rank 3'
714  read(*,*) casei
715 
716  select case(casei)
717 
718 
719 ! example 5 10
720 ! call of 5-point coefficients up to rank = 3
721  case (10)
722 
723  write(*,'(/(a))') ' Calculating 5-point coefficients'
724 
725  m02 = 1d2
726  m12 = 1d2
727  m22 = 3d2
728  m32 = 3d2
729  m42 = 4d2
730  masses2(0:4) = (/m02,m12,m22,m32,m42/)
731 
732  call getinvariants(5,mominv(1:10))
733  p10 = mominv(1)
734  p21 = mominv(2)
735  p32 = mominv(3)
736  p43 = mominv(4)
737  p40 = mominv(5)
738  p20 = mominv(6)
739  p31 = mominv(7)
740  p42 = mominv(8)
741  p30 = mominv(9)
742  p41 = mominv(10)
743 
744 ! allocate fields for tensor coefficients
745  rank = 3
746  allocate(ecoeff(0:rank/2,0:rank,0:rank,0:rank,0:rank))
747  allocate(ecoeffuv(0:rank/2,0:rank,0:rank,0:rank,0:rank))
748  allocate(eerr(0:rank))
749 
750 ! Ecoeff(n0,n1,n2,n3,n4) contains the results for the corresponding tensor coefficent
751 ! Ecoeffuv(n0,n1,n2,n3,n4) contains the results for the UV-singular part of the tensor coefficent
752 ! (coefficient of 2/(4-D)=1/eps term)
753 ! Eerr(r) is optional and contains an absolute estimate of the error of the
754 ! tensor coefficients of rank r (not involving metric tensors)
755  call e_cll(ecoeff,ecoeffuv,p10,p21,p32,p43,p40,p20,p31,p42,p30,p41,m02,m12,m22,m32,m42,rank,eerr)
756 
757 ! equivalent call
758  call e_cll(ecoeff,ecoeffuv,mominv(1:10),masses2(0:4),rank,eerr)
759 
760 ! equivalent calls, but all coefficients are in 1-dimensional array TNcoeff
761  n = 5
762  allocate(tncoeff(getnc_cll(n,rank)))
763  allocate(tncoeffuv(getnc_cll(n,rank)))
764  call e_cll(tncoeff,tncoeffuv,p10,p21,p32,p43,p40,p20,p31,p42,p30,p41,m02,m12,m22,m32,m42,rank,eerr)
765 ! or
766  call e_cll(tncoeff,tncoeffuv,mominv(1:10),masses2(0:4),rank,eerr)
767 ! or
768  call tn_cll(tncoeff,tncoeffuv,mominv(1:10),masses2(0:4),n,rank,eerr)
769 
770 ! writes the result to 'demo_5point_example10.dat'
771  call writeresulte(casen,casei,ecoeff,ecoeffuv,mominv(1:10),masses2(0:4),rank,eerr)
772 
773  deallocate(ecoeff,ecoeffuv,eerr,tncoeff,tncoeffuv)
774 
775 ! example 5 20
776 ! call of 5-point tensor components up to rank = 3
777  case (20)
778 
779  write(*,'(/(a))') ' Calculating 5-point tensor components'
780 
781  m02 = 1d2
782  m12 = 1d2
783  m22 = 3d2
784  m32 = 3d2
785  m42 = 4d2
786  masses2(0:4) = (/m02,m12,m22,m32,m42/)
787 
788  call getinvariants(5,mominv(1:10),momvec(0:3,1:4))
789  p1vec = momvec(0:3,1)
790  p2vec = momvec(0:3,2)
791  p3vec = momvec(0:3,3)
792  p4vec = momvec(0:3,4)
793  p10 = mominv(1)
794  p21 = mominv(2)
795  p32 = mominv(3)
796  p43 = mominv(4)
797  p40 = mominv(5)
798  p20 = mominv(6)
799  p31 = mominv(7)
800  p42 = mominv(8)
801  p30 = mominv(9)
802  p41 = mominv(10)
803 
804 ! allocate fields for tensor components
805  n = 5
806  rank = 3
807  allocate(eten(0:rank,0:rank,0:rank,0:rank))
808  allocate(etenuv(0:rank,0:rank,0:rank,0:rank))
809  allocate(eerr(0:rank))
810 
811 ! Eten(n0,n1,n2,n3) contains the results for the corresponding tensor component
812 ! Etenuv(n0,n1,n2,n3) contains the results for the UV-singular part of the tensor component
813 ! (coefficient of 2/(4-D)=1/eps term)
814 ! Eerr(r) is optional and contains an absolute estimate of the error of the
815 ! tensor coefficients of rank r (not involving metric tensors)
816  call eten_cll(eten,etenuv,p1vec,p2vec,p3vec,p4vec,p10,p21,p32,p43,p40,p20,p31,p42,p30,p41,m02,m12,m22,m32,m42,rank,eerr)
817 
818 ! equivalent calls
819  call eten_cll(eten,etenuv,momvec(0:3,1:4),mominv(1:10),masses2(0:4),rank,eerr)
820 ! or
821  call tnten_cll(eten,etenuv,momvec(0:3,1:4),mominv(1:10),masses2(0:4),n,rank,eerr)
822 
823 ! equivalent calls, but with all components are in 1-dimensional array TNten
824  allocate(tnten(getnt_cll(rank)))
825  allocate(tntenuv(getnt_cll(rank)))
826  call eten_cll(tnten,tntenuv,p1vec,p2vec,p3vec,p4vec,p10,p21,p32,p43,p40,p20,p31,p42,p30,p41,m02,m12,m22,m32,m42,rank,eerr)
827 ! or
828  call eten_cll(tnten,tntenuv,momvec(0:3,1:4),mominv(1:10),masses2(0:4),rank,eerr)
829 ! or
830  call tnten_cll(tnten,tntenuv,momvec(0:3,1:4),mominv(1:10),masses2(0:4),n,rank,eerr)
831 
832 ! writes the result to 'demo_5point_example20.dat'
833  call writeresulteten(casen,casei,eten,etenuv,momvec(0:3,1:4),mominv(1:10),masses2(0:4),rank,eerr)
834 
835  deallocate(eten,etenuv,eerr,tnten,tntenuv)
836 
837  case default
838 
839  write(*,'(/(a))') ' Input does not correspond to predefined sample'
840 
841  end select
842 
843 ! 6-point function calls
844  case(6)
845 
846  write(*,*)
847  write(*,*) 'Choose test case for 6-point function:'
848  write (*,*) '10) tensor 6-point function: coefficients up to rank 3'
849  write (*,*) '20) tensor 6-point function: tensor components of rank 3'
850  read(*,*) casei
851 
852  select case(casei)
853 
854 ! example 6 10
855 ! call of 6-point coefficients up to rank = 3
856  case (10)
857 
858  write(*,'(/(a))') ' Calculating 6-point coefficients'
859 
860  m02 = 1d2
861  m12 = 1d2
862  m22 = 3d2
863  m32 = 3d2
864  m42 = 4d2
865  m52 = 5d2
866  masses2(0:5) = (/m02,m12,m22,m32,m42,m52/)
867 
868  call getinvariants(6,mominv(1:15))
869 
870  p10 = mominv(1)
871  p21 = mominv(2)
872  p32 = mominv(3)
873  p43 = mominv(4)
874  p54 = mominv(5)
875  p50 = mominv(6)
876  p20 = mominv(7)
877  p31 = mominv(8)
878  p42 = mominv(9)
879  p53 = mominv(10)
880  p40 = mominv(11)
881  p51 = mominv(12)
882  p30 = mominv(13)
883  p41 = mominv(14)
884  p52 = mominv(15)
885 
886 ! allocate fields for tensor coefficients
887  rank = 3
888  allocate(fcoeff(0:rank/2,0:rank,0:rank,0:rank,0:rank,0:rank))
889  allocate(fcoeffuv(0:rank/2,0:rank,0:rank,0:rank,0:rank,0:rank))
890  allocate(ferr(0:rank))
891 
892 ! Fcoeff(n0,n1,n2,n3,n4) contains the results for the corresponding tensor coefficent
893 ! Fcoeffuv(n0,n1,n2,n3,n4) contains the results for the UV-singular part of the tensor coefficent
894 ! (coefficient of 2/(4-D)=1/eps term)
895 ! Ferr(r) is optional and contains an absolute estimate of the error of the
896 ! tensor coefficients of rank r (not involving metric tensors)
897  call f_cll(fcoeff,fcoeffuv,p10,p21,p32,p43,p54,p50,p20,p31,p42,p53,p40, &
898  p51,p30,p41,p52,m02,m12,m22,m32,m42,m52,rank,ferr)
899 
900 ! equivalent call
901  call f_cll(fcoeff,fcoeffuv,mominv(1:15),masses2(0:5),rank,ferr)
902 
903 ! equivalent calls, but all coefficients are in 1-dimensional array TNcoeff
904  n = 6
905  allocate(tncoeff(getnc_cll(n,rank)))
906  allocate(tncoeffuv(getnc_cll(n,rank)))
907  call f_cll(tncoeff,tncoeffuv,p10,p21,p32,p43,p54,p50,p20,p31,p42,p53,p40, &
908  p51,p30,p41,p52,m02,m12,m22,m32,m42,m52,rank,ferr)
909 ! or
910  call f_cll(tncoeff,tncoeffuv,mominv(1:15),masses2(0:5),rank,ferr)
911 ! or
912  call tn_cll(tncoeff,tncoeffuv,mominv(1:15),masses2(0:5),n,rank,ferr)
913 
914 
915 ! writes the result to 'demo_6point_example10.dat'
916  call writeresultf(casen,casei,fcoeff,fcoeffuv,mominv(1:15),masses2(0:5),rank,ferr)
917 
918  deallocate(fcoeff,fcoeffuv,ferr,tncoeff,tncoeffuv)
919 
920 ! example 6 20
921 ! call of 6-point tensor components up to rank = 3
922  case (20)
923 
924  write(*,'(/(a))') ' Calculating 6-point tensor components'
925 
926  m02 = 1d2
927  m12 = 1d2
928  m22 = 3d2
929  m32 = 3d2
930  m42 = 4d2
931  m52 = 5d2
932  masses2(0:5) = (/m02,m12,m22,m32,m42,m52/)
933 
934  call getinvariants(6,mominv(1:15),momvec(0:3,1:5))
935 
936  p1vec = momvec(0:3,1)
937  p2vec = momvec(0:3,2)
938  p3vec = momvec(0:3,3)
939  p4vec = momvec(0:3,4)
940  p5vec = momvec(0:3,5)
941  p10 = mominv(1)
942  p21 = mominv(2)
943  p32 = mominv(3)
944  p43 = mominv(4)
945  p54 = mominv(5)
946  p50 = mominv(6)
947  p20 = mominv(7)
948  p31 = mominv(8)
949  p42 = mominv(9)
950  p53 = mominv(10)
951  p40 = mominv(11)
952  p51 = mominv(12)
953  p30 = mominv(13)
954  p41 = mominv(14)
955  p52 = mominv(15)
956 
957 ! allocate fields for tensor coefficients
958  n = 6
959  rank = 3
960  allocate(ften(0:rank,0:rank,0:rank,0:rank))
961  allocate(ftenuv(0:rank,0:rank,0:rank,0:rank))
962  allocate(ferr(0:rank))
963 
964 ! Ften(n0,n1,n2,n3) contains the results for the corresponding tensor component
965 ! Ftenuv(n0,n1,n2,n3) contains the results for the UV-singular part of the tensor component
966 ! (coefficient of 2/(4-D)=1/eps term)
967 ! Ferr(r) is optional and contains an absolute estimate of the error of the
968 ! tensor coefficients of rank r (not involving metric tensors)
969  call ften_cll(ften,ftenuv,p1vec,p2vec,p3vec,p4vec,p5vec,p10,p21,p32,p43,p54,p50,p20,p31,p42,p53,p40, &
970  p51,p30,p41,p52,m02,m12,m22,m32,m42,m52,rank,ferr)
971 
972 ! equivalent calls
973  call ften_cll(ften,ftenuv,momvec(0:3,1:5),mominv(1:15),masses2(0:5),rank,ferr)
974 ! or
975  call tnten_cll(ften,ftenuv,momvec(0:3,1:5),mominv(1:15),masses2(0:5),n,rank,ferr)
976 
977 ! equivalent calls, but with all components are in 1-dimensional array TNten
978  allocate(tnten(getnt_cll(rank)))
979  allocate(tntenuv(getnt_cll(rank)))
980  call ften_cll(tnten,tntenuv,p1vec,p2vec,p3vec,p4vec,p5vec,p10,p21,p32,p43,p54,p50,p20,p31,p42,p53,p40, &
981  p51,p30,p41,p52,m02,m12,m22,m32,m42,m52,rank,ferr)
982 ! or
983  call ften_cll(tnten,tntenuv,momvec(0:3,1:5),mominv(1:15),masses2(0:5),rank,ferr)
984 ! or
985  call tnten_cll(tnten,tntenuv,momvec(0:3,1:5),mominv(1:15),masses2(0:5),n,rank,ferr)
986 
987 ! writes the result to 'demo_6point_example10.dat'
988  call writeresultften(casen,casei,ften,ftenuv,momvec(0:3,1:5),mominv(1:15),masses2(0:5),rank,ferr)
989 
990  deallocate(ften,ftenuv,ferr,tnten,tntenuv)
991 
992  case default
993 
994  write(*,'(/(a))') ' Input does not correspond to predefined sample'
995 
996  end select
997 
998 ! calls for derivatives of scalar 2-point function
999  case(12)
1000 
1001  write(*,*)
1002  write(*,*) 'Choose test case for 2-point function derivative:'
1003  write (*,*) ' 1) scalar 2-point function derivative'
1004  write (*,*) ' 2) vector 2-point function derivative'
1005  write (*,*) ' 3) tensor 2-point function derivative'
1006  write (*,*) ' 4) mass-singular scalar 2-point function derivative'
1007  write (*,*) '10) tensor 2-point function derivatives up to rank 2'
1008  read(*,*) casei
1009 
1010  select case(casei)
1011 
1012 ! example 12 1
1013 ! call of derivative of scalar 2-point function
1014  case(1)
1015 
1016  write(*,'(/(a)/)') ' Calculating derivative of scalar 2-point function'
1017 
1018  p10 = 1d3
1019  m02 = 4d2
1020  m12 = 1d2
1021 
1022  mominv(1:1) = p10
1023  masses2(0:1) = (/m02,m12/)
1024 
1025 ! DB0 contains the results for the derivative of the scalar 2-point function
1026  call db0_cll(db0,p10,m02,m12)
1027 
1028 ! equivalent call
1029  call db0_cll(db0,mominv(1:1),masses2(0:1))
1030 
1031 ! writes the result to 'demo_2point_derivative_example01.dat'
1032  call writeresultdb0(casen,casei,db0,mominv(1:1),masses2(0:1))
1033 
1034 
1035 ! example 12 2
1036 ! call of derivative of vector 2-point function coefficient
1037  case(2)
1038 
1039  write(*,'(/(a)/)') ' Calculating derivative of vector 2-point function coefficient'
1040 
1041  p10 = 1d3
1042  m02 = 4d2
1043  m12 = 1d2
1044 
1045  mominv(1:1) = p10
1046  masses2(0:1) = (/m02,m12/)
1047 
1048 ! DB1 contains the results for the derivative of the scalar 2-point function
1049  call db1_cll(db1,p10,m02,m12)
1050 
1051 ! equivalent call
1052  call db1_cll(db1,mominv(1:1),masses2(0:1))
1053 
1054 ! writes the result to 'demo_2point_derivative_example02.dat'
1055  rank = 1
1056  allocate(dbcoeff(0:rank/2,0:rank))
1057  allocate(dbcoeffuv(0:rank/2,0:rank))
1058  allocate(dberr(0:rank))
1059  dbcoeff(0,1) = db1
1060  call writeresultdb(casen,casei,dbcoeff,dbcoeffuv,mominv(1:1),masses2(0:1),1,1)
1061 
1062 
1063 ! example 12 3
1064 ! call of derivative of tensor 2-point function coefficient
1065  case(3)
1066 
1067  write(*,'(/(a)/)') ' Calculating derivative of tensor 2-point function coefficient'
1068 
1069  p10 = 1d3
1070  m02 = 4d2
1071  m12 = 1d2
1072 
1073  mominv(1:1) = p10
1074  masses2(0:1) = (/m02,m12/)
1075 
1076 ! DB00 contains the results for the derivative of the scalar 2-point function
1077  call db00_cll(db00,db00uv,p10,m02,m12)
1078 
1079 ! equivalent call
1080  call db00_cll(db00,db00uv,mominv(1:1),masses2(0:1))
1081 
1082 ! writes the result to 'demo_2point_derivative_example03.dat'
1083  rank = 2
1084  allocate(dbcoeff(0:rank/2,0:rank))
1085  allocate(dbcoeffuv(0:rank/2,0:rank))
1086  allocate(dberr(0:rank))
1087  dbcoeff(1,0) = db00
1088  call writeresultdb(casen,casei,dbcoeff,dbcoeffuv,mominv(1:1),masses2(0:1),2,2)
1089 
1090 
1091 ! example 12 4
1092 ! call of derivative of mass-singular scalar 2-point function
1093  case(4)
1094 
1095  write(*,'(/(a))') ' Calculating derivative of mass-singular scalar 2-point function'
1096  write(*,'((a)/)') ' in mass regularization'
1097  write(*,'((a)/)') ' mass-singularity described by log(m12)'
1098 
1099  p10 = 0d0
1100  m02 = 1d-4
1101  m12 = 1d-4
1102 
1103 ! clears list of regulator masses
1104  call clearminf2_cll
1105 ! adds m12 to list of regulator masses
1106 ! these are treated as infinitesimal apart from their appearance
1107 ! in IR-singular logarithms, where the actual values are taken
1108  call addminf2_cll(m12)
1109 
1110  mominv(1:1) = p10
1111  masses2(0:1) = (/m02,m12/)
1112 
1113 ! DB0 contains the results for the scalar 2-point function
1114  call db0_cll(db0,p10,m02,m12)
1115 
1116 ! equivalent call
1117  call db0_cll(db0,mominv(1:1),masses2(0:1))
1118 
1119 ! writes the result to 'demo_2point_example04.dat'
1120  call writeresultdb0(casen,casei,db0,mominv(1:1),masses2(0:1))
1121 
1122 
1123 ! example 12 10
1124 ! call of derivatives of tensor 2-point function coefficients up to rank = 2
1125  case (10)
1126 
1127  write(*,'(/(a))') ' Calculating derivatives of tensor 2-point coefficients'
1128 
1129  p10 = 1d3
1130  m02 = 4d2
1131  m12 = 1d2
1132 
1133  mominv(1:1) = (/p10/)
1134  masses2(0:1) = (/m02,m12/)
1135 
1136 ! allocate fields for tensor coefficients
1137  rank = 2
1138  allocate(dbcoeff(0:rank/2,0:rank))
1139  allocate(dbcoeffuv(0:rank/2,0:rank))
1140  allocate(dberr(0:rank))
1141 
1142 
1143 ! DBcoeff(n0,n1) contains the results for the corresponding tensor coefficent derivative
1144 ! DBcoeffuv(n0,n1) contains the results for the UV-singular part of the tensor coefficent derivative
1145 ! (coefficient of 2/(4-D)=1/eps term)
1146 ! DBerr(r) is optional and contains an absolute estimate of the error of the derivative of the
1147 ! tensor coefficients of rank r (not involving metric tensors)
1148  call db_cll(dbcoeff,dbcoeffuv,p10,m02,m12,rank,dberr)
1149 
1150 ! equivalent call
1151  call db_cll(dbcoeff,dbcoeffuv,mominv(1:1),masses2(0:1),rank,dberr)
1152 
1153 ! writes the result to 'demo_2point_derivative_example10.dat'
1154  call writeresultdb(casen,casei,dbcoeff,dbcoeffuv,mominv(1:1),masses2(0:1),0,rank,berr)
1155 
1156  deallocate(dbcoeff,dbcoeffuv,dberr)
1157 
1158  case default
1159 
1160  write(*,'(/(a))') ' Input does not correspond to predefined sample'
1161 
1162  end select
1163 
1164  case default
1165 
1166  write(*,'(/(a))') ' Input does not correspond to predefined sample'
1167 
1168  end select
1169 
1170  contains
1171 
1172 
1173  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1174  ! subroutine writeresultA(caseN,casei,Acoeff,Acoeffuv,masses2,rank,Aerr)
1175  !
1176  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1177 
1178  subroutine writeresulta(caseN,casei,Acoeff,Acoeffuv,masses2,rank,Aerr)
1179  implicit none
1180 
1181  integer, intent(in) :: rank,caseN,casei
1182  double complex, intent(in) :: masses2(0:0)
1183  double complex, intent(in) :: Acoeff(0:rank/2)
1184  double complex, intent(in) :: Acoeffuv(0:rank/2)
1185  double precision, optional, intent(in) :: Aerr(0:rank)
1186  integer :: r,n0,mode
1187  character(len=99) :: fname
1188  character(len=*),parameter :: fmt1 = "(A9,' = ',es23.16,' + i*',es23.16)"
1189  character(len=*),parameter :: fmt2 = "(A9,' = ',es23.16)"
1190  character(len=*),parameter :: fmt10 = "(' Acoeff(',i1,') = ',es23.16,' + i*',es23.16)"
1191  character(len=*),parameter :: fmt11 = "(' Aerr(',i1,') = ',es23.16)"
1192  character(len=*),parameter :: fmt12 = "(' A0 = ',es23.16,' + i*',es23.16)"
1193  character(len=*),parameter :: fmt13 = "(' Acoeffuv(',i1,') = ',es23.16,' + i*',es23.16)"
1194 
1195 ! output file for result
1196  call getmode_cll(mode)
1197  select case (mode)
1198  case(1)
1199  fname = 'demo_1point_example00_coli.dat'
1200  case(2)
1201  fname = 'demo_1point_example00_dd.dat'
1202  case(3)
1203  fname = 'demo_1point_example00_comp.dat'
1204  end select
1205  write(fname(20:21),'(i2.2)') casei
1206 
1207  open(unit=50,file=trim(fname),status='unknown')
1208 
1209  call getmuuv2_cll(muuv2)
1210 ! call GetMuIR2_cll(MuIR2)
1211  call getdeltauv_cll(deltauv)
1212 ! call GetdeltaIR_cll(DeltaIR1,DeltaIR2)
1213 
1214  write (50,'(a37,i2,i3/)') ' Result for 1-point function, example',casen,casei
1215  write (50,'((a))') ' '
1216  write (50,'((a))') ' ------- '
1217  write (50,'((a))') ' / \ '
1218  write (50,'((a))') ' / \ '
1219  write (50,'((a))') ' 0d0 ----- 0 | m02 '
1220  write (50,'((a))') ' \ / '
1221  write (50,'((a))') ' \ / '
1222  write (50,'((a))') ' ------- '
1223  write (50,'((a))') ' '
1224  write (50,'((a))') ''
1225  write (50,'((a))') ' Input:'
1226  write (50,fmt1) ' m02 ',masses2(0)
1227  write (50,fmt2) ' muUV2 ',muuv2
1228 ! write (50,fmt2) ' muIR2 ',muIR2
1229  write (50,fmt2) ' deltaUV ',deltauv
1230 ! write (50,fmt2) ' deltaIR1',deltaIR2
1231 ! write (50,fmt2) ' deltaIR2',deltaIR1
1232  write (50,'((a))') ''
1233  write (50,'((a))') ' Conventions:'
1234  write (50,'((a))') ''
1235  if(rank.gt.0) then
1236  write (50,'((a))') ' (2*pi*mu)^(4-D) '
1237  write (50,'((a))') ' A = --------------- \int d^D q f(q)'
1238  write (50,'((a))') ' i*pi^2 '
1239 
1240  write (50,'((a))')
1241  write (50,'((a))') ' = A_fin(muUV2) + a_UV*DeltaUV '
1242  else
1243  write (50,'((a))') ' (2*pi*mu)^(4-D) '
1244  write (50,'((a))') ' A0 = --------------- \int d^D q f(q)'
1245  write (50,'((a))') ' i*pi^2 '
1246 
1247  write (50,'((a))')
1248  write (50,'((a))') ' = A0_fin(muUV2) + a_UV*DeltaUV '
1249  end if
1250  write (50,'((a))')
1251  write (50,'((a))') ' where'
1252  write (50,'((a))')
1253  write (50,'((a))') &
1254  ' c(epsUV) '
1255  write (50,'((a))') &
1256  ' DeltaUV = -------- '
1257  write (50,'((a))') &
1258  ' epsUV '
1259  write (50,'((a))')
1260  write (50,'((a))') ' c(eps) = (4*pi)^eps\Gamma(1+eps), D = 4 -2*eps '
1261  write (50,'((a))')
1262  write (50,'((a))') ' you can freely choose the regularization parameters'
1263  write (50,'((a))') ' of UV origin: muUV2 = mu^2, DeltaUV '
1264  write (50,'((a))')
1265  write (50,'((a))') ' note:'
1266  write (50,'((a))') ' - we effectively factor out a factor c(eps) '
1267  write (50,'((a))') ' - by default DeltaUV = 0 '
1268 
1269  write (50,'((a))') ''
1270  write (50,'((a))') ' Results:'
1271 
1272  if (rank.gt.0) then
1273  do n0=0,rank/2
1274  write (50,fmt10) n0,acoeff(n0)
1275  end do
1276 
1277  write (50,'(/(a))') ' Error estimates:'
1278  do r=0,rank
1279  write (50,fmt11) r,aerr(r)
1280  end do
1281  else
1282  write (50,fmt12) acoeff(0)
1283  end if
1284 
1285  write(*,'(/(a),(a)/)') ' The result has been written to the file ' &
1286  ,trim(fname)
1287 
1288  end subroutine writeresulta
1289 
1290 
1291  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1292  ! subroutine writeresultA0(caseN,casei,A0,masses2)
1293  !
1294  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1295 
1296  subroutine writeresulta0(caseN,casei,A0,masses2)
1297  implicit none
1298 
1299  integer, intent(in) :: caseN,casei
1300  double complex, intent(in) :: A0
1301  double complex, intent(in) :: masses2(0:0)
1302  double complex :: Acoeff(0:0)
1303  double complex :: Acoeffuv(0:0)
1304 
1305 
1306  acoeff(0:0) = a0
1307  acoeffuv(0:0) = masses2
1308  call writeresulta(casen,casei,acoeff,acoeffuv,masses2,0)
1309 
1310  end subroutine writeresulta0
1311 
1312 
1313  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1314  ! subroutine writeresultB(caseN,casei,Bcoeff,Bcoeffuv,MomInv,masses2,rank,Berr)
1315  !
1316  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1317 
1318  subroutine writeresultb(caseN,casei,Bcoeff,Bcoeffuv,MomInv,masses2,rank,Berr)
1319  implicit none
1320 
1321  integer, intent(in) :: rank,caseN,casei
1322  double complex, intent(in) :: MomInv(1), masses2(0:1)
1323  double complex, intent(in) :: Bcoeff(0:rank/2,0:rank)
1324  double complex, intent(in) :: Bcoeffuv(0:rank/2,0:rank)
1325  double precision, optional, intent(in) :: Berr(0:rank)
1326  integer :: r,n0,n1,mode
1327  integer,parameter :: rankuv=0
1328  character(len=99) :: fname
1329  character(len=*),parameter :: fmt1 = "(A9,' = ',es23.16,' + i*',es23.16)"
1330  character(len=*),parameter :: fmt2 = "(A9,' = ',es23.16)"
1331  character(len=*),parameter :: fmt10 = "(' Bcoeff(',i1,1(',',i1),') = ',es23.16,' + i*',es23.16)"
1332  character(len=*),parameter :: fmt11 = "(' Berr(',i1,') = ',es23.16)"
1333  character(len=*),parameter :: fmt12 = "(' B0 = ',es23.16,' + i*',es23.16)"
1334  character(len=*),parameter :: fmt13 = "(' Bcoeffuv(',i1,1(',',i1),') = ',es23.16,' + i*',es23.16)"
1335 
1336 ! output file for result
1337  call getmode_cll(mode)
1338  select case (mode)
1339  case(1)
1340  fname = 'demo_2point_example00_coli.dat'
1341  case(2)
1342  fname = 'demo_2point_example00_dd.dat'
1343  case(3)
1344  fname = 'demo_2point_example00_comp.dat'
1345  end select
1346  write(fname(20:21),'(i2.2)') casei
1347 
1348  open(unit=50,file=trim(fname),status='unknown')
1349 
1350  call getmuuv2_cll(muuv2)
1351  call getmuir2_cll(muir2)
1352  call getdeltauv_cll(deltauv)
1353  call getdeltair_cll(deltair1,deltair2)
1354 
1355  write (50,'(a37,i2,i3/)') ' Result for 2-point function, example',casen,casei
1356  write (50,'((a))') ' m12 '
1357  write (50,'((a))') ' ------- '
1358  write (50,'((a))') ' / 1 \ '
1359  write (50,'((a))') ' / \ '
1360  write (50,'((a))') ' p10 ----- ----- p10 '
1361  write (50,'((a))') ' \ / '
1362  write (50,'((a))') ' \ 0 / '
1363  write (50,'((a))') ' ------- '
1364  write (50,'((a))') ' m02 '
1365  write (50,'((a))') ''
1366  write (50,'((a))') ' Input:'
1367  write (50,fmt1) ' p10 ',mominv(1)
1368  write (50,fmt1) ' m02 ',masses2(0)
1369  write (50,fmt1) ' m12 ',masses2(1)
1370  if (rank.ge.rankuv) then
1371  write (50,fmt2) ' muUV2 ',muuv2
1372  end if
1373 ! write (50,fmt2) ' muIR2 ',muIR2
1374  if (rank.ge.rankuv) then
1375  write (50,fmt2) ' deltaUV ',deltauv
1376  end if
1377 ! write (50,fmt2) ' deltaIR1',deltaIR2
1378 ! write (50,fmt2) ' deltaIR2',deltaIR1
1379  write (50,'((a))') ''
1380  write (50,'((a))') ' Conventions:'
1381  write (50,'((a))') ''
1382  if(rank.gt.0) then
1383  write (50,'((a))') ' (2*pi*mu)^(4-D) '
1384  write (50,'((a))') ' B = --------------- \int d^D q f(q,p_i)'
1385  write (50,'((a))') ' i*pi^2 '
1386 
1387  write (50,'((a))')
1388  write (50,'((a))') ' = B_fin(muUV2,muIR2) + a_UV*DeltaUV '
1389  else
1390  write (50,'((a))') ' (2*pi*mu)^(4-D) '
1391  write (50,'((a))') ' B0 = --------------- \int d^D q f(q,p_i)'
1392  write (50,'((a))') ' i*pi^2 '
1393 
1394  write (50,'((a))')
1395  write (50,'((a))') ' = B0_fin(muUV2,muIR2) + a_UV*DeltaUV '
1396  end if
1397  write (50,'((a))')
1398  write (50,'((a))') ' where'
1399  write (50,'((a))')
1400  write (50,'((a))') &
1401  ' c(epsUV) '
1402  write (50,'((a))') &
1403  ' DeltaUV = -------- '
1404  write (50,'((a))') &
1405  ' epsUV '
1406  write (50,'((a))')
1407  write (50,'((a))') ' c(eps) = (4*pi)^eps\Gamma(1+eps), D = 4 -2*eps '
1408  write (50,'((a))')
1409  write (50,'((a))') ' you can freely choose the regularization parameters'
1410  write (50,'((a))') ' of UV origin: muUV2 = mu^2, DeltaUV '
1411  write (50,'((a))')
1412  write (50,'((a))') ' note:'
1413  write (50,'((a))') ' - we effectively factor out a factor c(eps) '
1414  write (50,'((a))') ' - by default DeltaUV = 0 '
1415 
1416  write (50,'((a))') ''
1417  write (50,'((a))') ' Results:'
1418 
1419  if (rank.gt.0) then
1420  do r=0,rank
1421  do n0=0,r/2
1422  n1=r-2*n0
1423  write (50,fmt10) n0,n1,bcoeff(n0,n1)
1424  end do
1425  end do
1426 
1427 ! do r=0,rank
1428 ! do n0=0,r/2
1429 ! n1=r-2*n0
1430 ! write (50,fmt13) n0,n1,Bcoeffuv(n0,n1)
1431 ! end do
1432 ! end do
1433 !
1434  write (50,'(/(a))') ' Error estimates:'
1435  do r=0,rank
1436  write (50,fmt11) r,berr(r)
1437  end do
1438  else
1439  write (50,fmt12) bcoeff(0,0)
1440  end if
1441 
1442  write(*,'(/(a),(a)/)') ' The result has been written to the file ' &
1443  ,trim(fname)
1444 
1445  end subroutine writeresultb
1446 
1447 
1448  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1449  ! subroutine writeresultB0(caseN,casei,B0,MomInv,masses2)
1450  !
1451  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1452 
1453  subroutine writeresultb0(caseN,casei,B0,MomInv,masses2)
1454  implicit none
1455 
1456  integer, intent(in) :: caseN,casei
1457  double complex, intent(in) :: B0
1458  double complex, intent(in) :: MomInv(1), masses2(0:1)
1459  double complex :: Bcoeff(0:0,0:0)
1460  double complex :: Bcoeffuv(0:0,0:0)
1461 
1462 
1463  bcoeff(0:0,0:0) = b0
1464  bcoeffuv(0:0,0:0) = 1d0
1465  call writeresultb(casen,casei,bcoeff,bcoeffuv,mominv,masses2,0)
1466 
1467  end subroutine writeresultb0
1468 
1469 
1470 
1471  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1472  ! subroutine writeresultC(caseN,casei,Ccoeff,Ccoeffuv,MomInv,masses2,rank,Cerr)
1473  !
1474  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1475 
1476  subroutine writeresultc(caseN,casei,Ccoeff,Ccoeffuv,MomInv,masses2,rank,Cerr)
1477  implicit none
1478 
1479  integer, intent(in) :: rank,caseN,casei
1480  double complex, intent(in) :: MomInv(3), masses2(0:2)
1481  double complex, intent(in) :: Ccoeff(0:rank/2,0:rank,0:rank)
1482  double complex, intent(in) :: Ccoeffuv(0:rank/2,0:rank,0:rank)
1483  double precision, optional, intent(in) :: Cerr(0:rank)
1484  integer :: r,n0,n1,n2,mode
1485  integer,parameter :: rankuv=2
1486  character(len=99) :: fname
1487  character(len=*),parameter :: fmt1 = "(A9,' = ',es23.16,' + i*',es23.16)"
1488  character(len=*),parameter :: fmt2 = "(A9,' = ',es23.16)"
1489  character(len=*),parameter :: fmt10 = "(' Ccoeff(',i1,2(',',i1),') = ',es23.16,' + i*',es23.16)"
1490  character(len=*),parameter :: fmt11 = "(' Cerr(',i1,') = ',es23.16)"
1491  character(len=*),parameter :: fmt12 = "(' C0 = ',es23.16,' + i*',es23.16)"
1492 
1493 ! output file for result
1494  call getmode_cll(mode)
1495  select case (mode)
1496  case(1)
1497  fname = 'demo_3point_example00_coli.dat'
1498  case(2)
1499  fname = 'demo_3point_example00_dd.dat'
1500  case(3)
1501  fname = 'demo_3point_example00_comp.dat'
1502  end select
1503  write(fname(20:21),'(i2.2)') casei
1504 
1505 ! write(*,*) casei,fname
1506 
1507  open(unit=50,file=trim(fname),status='unknown')
1508 
1509  call getmuuv2_cll(muuv2)
1510  call getmuir2_cll(muir2)
1511  call getdeltauv_cll(deltauv)
1512  call getdeltair_cll(deltair1,deltair2)
1513 
1514  write (50,'(a37,i2,i3/)') ' Result for 3-point function, example',casen,casei
1515  write (50,'(a63,i2,i3,a)') ' The corresponding code can be found in demo.f90 under ''example',casen,casei,''''
1516  write (50,'((a))')
1517  write (50,'((a))') ' p21 '
1518  write (50,'((a))') ' | '
1519  write (50,'((a))') ' | '
1520  write (50,'((a))') ' / \ '
1521  write (50,'((a))') ' / \ '
1522  write (50,'((a))') ' m12 /1 2\ m22 '
1523  write (50,'((a))') ' / \ '
1524  write (50,'((a))') ' / 0 \ '
1525  write (50,'((a))') ' p10 --------------------- p20 '
1526  write (50,'((a))') ' m02 '
1527  write (50,'((a))') ''
1528  write (50,'((a))') ' Input:'
1529  write (50,fmt1) ' p10 ',mominv(1)
1530  write (50,fmt1) ' p21 ',mominv(2)
1531  write (50,fmt1) ' p20 ',mominv(3)
1532  write (50,fmt1) ' m02 ',masses2(0)
1533  write (50,fmt1) ' m12 ',masses2(1)
1534  write (50,fmt1) ' m22 ',masses2(2)
1535  if (rank.ge.2) then
1536  write (50,fmt2) ' muUV2 ',muuv2
1537  end if
1538  write (50,fmt2) ' muIR2 ',muir2
1539  if (rank.ge.2) then
1540  write (50,fmt2) ' DeltaUV ',deltauv
1541  end if
1542  write (50,fmt2) ' DeltaIR1',deltair2
1543  write (50,fmt2) ' DeltaIR2',deltair1
1544  write (50,'((a))') ''
1545  write (50,'((a))') ' Conventions:'
1546  write (50,'((a))') ''
1547  if(rank.eq.0) then
1548  write (50,'((a))') ' (2*pi*mu)^(4-D) '
1549  write (50,'((a))') ' C0 = --------------- \int d^D q f(q,p_i)'
1550  write (50,'((a))') ' i*pi^2 '
1551 
1552  write (50,'((a))')
1553  write (50,'((a))') ' = C0_fin(muUV2,muIR2) + '// &
1554  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
1555  else if(rank.ge.rankuv) then
1556  write (50,'((a))') ' (2*pi*mu)^(4-D) '
1557  write (50,'((a))') ' C = --------------- \int d^D q f(q,p_i)'
1558  write (50,'((a))') ' i*pi^2 '
1559 
1560  write (50,'((a))')
1561  write (50,'((a))') ' = C_fin(muUV2,muIR2) + a_UV*DeltaUV + '// &
1562  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
1563  else
1564  write (50,'((a))') ' (2*pi*mu)^(4-D) '
1565  write (50,'((a))') ' C = --------------- \int d^D q f(q,p_i)'
1566  write (50,'((a))') ' i*pi^2 '
1567 
1568  write (50,'((a))')
1569  write (50,'((a))') ' = C_fin(muUV2,muIR2) + '// &
1570  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
1571  end if
1572  write (50,'((a))')
1573  write (50,'((a))') ' where'
1574  write (50,'((a))')
1575  if (rank.gt.0) then
1576  write (50,'((a))') &
1577  ' c(epsUV) c(epsIR) c(epsIR)'
1578  write (50,'((a))') &
1579  ' DeltaUV = --------, DeltaIR1 = --------, DeltaIR2 = --------'
1580  write (50,'((a))') &
1581  ' epsUV epsIR epsIR^2'
1582  else
1583  write (50,'((a))') &
1584  ' c(epsIR) c(epsIR)'
1585  write (50,'((a))') &
1586  ' DeltaIR1 = --------, DeltaIR2 = --------'
1587  write (50,'((a))') &
1588  ' epsIR epsIR^2'
1589  end if
1590  write (50,'((a))')
1591  write (50,'((a))') ' c(eps) = (4*pi)^eps\Gamma(1+eps), D = 4 -2*eps '
1592  write (50,'((a))')
1593  write (50,'((a))') ' you can freely choose the regularization parameters'
1594  if (rank.gt.0) then
1595  write (50,'((a))') ' of UV origin: muUV2 = mu^2, DeltaUV '
1596  end if
1597  write (50,'((a))') ' of IR origin: muIR2 = mu^2, DeltaIR1, DeltaIR2'
1598  write (50,'((a))')
1599  write (50,'((a))') ' note:'
1600  write (50,'((a))') ' - we effectively factor out a factor c(eps) '
1601  if (rank.gt.0) then
1602  write (50,'((a))') ' - by default DeltaUV = DeltaIR1 = DeltaIR2 = 0 '
1603  else
1604  write (50,'((a))') ' - by default DeltaIR1 = DeltaIR2 = 0 '
1605  end if
1606  write (50,'((a))') ' - suitable DeltaIR2 can be used to adapt the effective normalization'
1607 
1608  write (50,'((a))') ''
1609  write (50,'((a))') ' Results:'
1610 
1611  if (rank.gt.0) then
1612  do r=0,rank
1613  do n0=0,r/2
1614  do n1=0,r-2*n0
1615  n2 = r-2*n0-n1
1616  write (50,fmt10) n0,n1,n2,ccoeff(n0,n1,n2)
1617  end do
1618  end do
1619  end do
1620 
1621  write (50,'(/(a))') ' Error estimates:'
1622  do r=0,rank
1623  write (50,fmt11) r,cerr(r)
1624  end do
1625  else
1626  write (50,fmt12) ccoeff(0,0,0)
1627  end if
1628 
1629  write(*,'(/(a),(a)/)') ' The result has been written to the file ' &
1630  ,trim(fname)
1631 
1632  end subroutine writeresultc
1633 
1634 
1635  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1636  ! subroutine writeresultC0(caseN,casei,Ccoeff,Ccoeffuv,MomInv,masses2)
1637  !
1638  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1639 
1640  subroutine writeresultc0(caseN,casei,C0,MomInv,masses2)
1641  implicit none
1642 
1643  integer, intent(in) :: caseN,casei
1644  double complex, intent(in) :: C0
1645  double complex, intent(in) :: MomInv(3), masses2(0:2)
1646  double complex :: Ccoeff(0:0,0:0,0:0)
1647  double complex :: Ccoeffuv(0:0,0:0,0:0)
1648 
1649 
1650  ccoeff(0:0,0:0,0:0) = c0
1651  ccoeffuv(0:0,0:0,0:0) = 0d0
1652  call writeresultc(casen,casei,ccoeff,ccoeffuv,mominv,masses2,0)
1653 
1654  end subroutine writeresultc0
1655 
1656 
1657 
1658  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1659  ! subroutine writeresultCten(caseN,casei,Cten,Ctenuv,MomInv,masses2,rank,Cerr)
1660  !
1661  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1662 
1663  subroutine writeresultcten(caseN,casei,Cten,Ctenuv,MomVec,MomInv,masses2,rank,Cerr)
1664  implicit none
1665 
1666  integer, intent(in) :: rank,caseN,casei
1667  double complex, intent(in) :: MomVec(0:3,1:2), MomInv(3), masses2(0:2)
1668  double complex, intent(in) :: Cten(0:rank,0:rank,0:rank,0:rank)
1669  double complex, intent(in) :: Ctenuv(0:rank,0:rank,0:rank,0:rank)
1670  double precision, optional, intent(in) :: Cerr(0:rank)
1671  integer :: r,n0,n1,n2,n3,mode
1672  integer,parameter :: rankuv=2
1673  character(len=99) :: fname
1674  character(len=*),parameter :: fmt1 = "(A17,' = ',es23.16,' + i*',es23.16)"
1675  character(len=*),parameter :: fmt2 = "(A17,' = ',es23.16)"
1676  character(len=*),parameter :: fmt10 = "(' Cten(',i1,3(',',i1),') = ',es23.16,' + i*',es23.16)"
1677  character(len=*),parameter :: fmt11 = "(' Cerr(',i1,') = ',es23.16)"
1678  character(len=*),parameter :: fmt12 = "(' C0 = ',es23.16,' + i*',es23.16)"
1679 
1680 ! output file for result
1681  call getmode_cll(mode)
1682  select case (mode)
1683  case(1)
1684  fname = 'demo_3point_example00_coli.dat'
1685  case(2)
1686  fname = 'demo_3point_example00_dd.dat'
1687  case(3)
1688  fname = 'demo_3point_example00_comp.dat'
1689  end select
1690  write(fname(20:21),'(i2.2)') casei
1691 
1692 ! write(*,*) casei,fname
1693 
1694  open(unit=50,file=trim(fname),status='unknown')
1695 
1696  call getmuuv2_cll(muuv2)
1697  call getmuir2_cll(muir2)
1698  call getdeltauv_cll(deltauv)
1699  call getdeltair_cll(deltair1,deltair2)
1700 
1701  write (50,'(a37,i2,i3/)') ' Result for 3-point function, example',casen,casei
1702  write (50,'(a63,i2,i3,a)') ' The corresponding code can be found in demo.f90 under ''example',casen,casei,''''
1703  write (50,'((a))')
1704  write (50,'((a))') ' p21 '
1705  write (50,'((a))') ' | '
1706  write (50,'((a))') ' | '
1707  write (50,'((a))') ' / \ '
1708  write (50,'((a))') ' / \ '
1709  write (50,'((a))') ' m12,p1vec /1 2\ m22,p2vec '
1710  write (50,'((a))') ' / \ '
1711  write (50,'((a))') ' / 0 \ '
1712  write (50,'((a))') ' p10 --------------------- p20 '
1713  write (50,'((a))') ' m02 '
1714  write (50,'((a))') ''
1715  write (50,'((a))') ' Input:'
1716  write (50,fmt1) ' p1vec(0) ',momvec(0,1)
1717  write (50,fmt1) ' p1vec(1) ',momvec(1,1)
1718  write (50,fmt1) ' p1vec(2) ',momvec(2,1)
1719  write (50,fmt1) ' p1vec(3) ',momvec(3,1)
1720  write (50,fmt1) ' p2vec(0) ',momvec(0,2)
1721  write (50,fmt1) ' p2vec(1) ',momvec(1,2)
1722  write (50,fmt1) ' p2vec(2) ',momvec(2,2)
1723  write (50,fmt1) ' p2vec(3) ',momvec(3,2)
1724  write (50,fmt1) ' p10 ',mominv(1)
1725  write (50,fmt1) ' p21 ',mominv(2)
1726  write (50,fmt1) ' p20 ',mominv(3)
1727  write (50,fmt1) ' m02 ',masses2(0)
1728  write (50,fmt1) ' m12 ',masses2(1)
1729  write (50,fmt1) ' m22 ',masses2(2)
1730  if (rank.ge.2) then
1731  write (50,fmt2) ' muUV2 ',muuv2
1732  end if
1733  write (50,fmt2) ' muIR2 ',muir2
1734  if (rank.ge.2) then
1735  write (50,fmt2) ' DeltaUV ',deltauv
1736  end if
1737  write (50,fmt2) ' DeltaIR1 ',deltair2
1738  write (50,fmt2) ' DeltaIR2 ',deltair1
1739  write (50,'((a))') ''
1740  write (50,'((a))') ' Conventions:'
1741  write (50,'((a))') ''
1742  if(rank.eq.0) then
1743  write (50,'((a))') ' (2*pi*mu)^(4-D) '
1744  write (50,'((a))') ' C0 = --------------- \int d^D q f(q,p_i)'
1745  write (50,'((a))') ' i*pi^2 '
1746 
1747  write (50,'((a))')
1748  write (50,'((a))') ' = C0_fin(muUV2,muIR2) + '// &
1749  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
1750  else if(rank.ge.rankuv) then
1751  write (50,'((a))') ' (2*pi*mu)^(4-D) '
1752  write (50,'((a))') ' C = --------------- \int d^D q f(q,p_i)'
1753  write (50,'((a))') ' i*pi^2 '
1754 
1755  write (50,'((a))')
1756  write (50,'((a))') ' = C_fin(muUV2,muIR2) + a_UV*DeltaUV + '// &
1757  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
1758  else
1759  write (50,'((a))') ' (2*pi*mu)^(4-D) '
1760  write (50,'((a))') ' C = --------------- \int d^D q f(q,p_i)'
1761  write (50,'((a))') ' i*pi^2 '
1762 
1763  write (50,'((a))')
1764  write (50,'((a))') ' = C_fin(muUV2,muIR2) + '// &
1765  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
1766  end if
1767  write (50,'((a))')
1768  write (50,'((a))') ' where'
1769  write (50,'((a))')
1770  if (rank.gt.0) then
1771  write (50,'((a))') &
1772  ' c(epsUV) c(epsIR) c(epsIR)'
1773  write (50,'((a))') &
1774  ' DeltaUV = --------, DeltaIR1 = --------, DeltaIR2 = --------'
1775  write (50,'((a))') &
1776  ' epsUV epsIR epsIR^2'
1777  else
1778  write (50,'((a))') &
1779  ' c(epsIR) c(epsIR)'
1780  write (50,'((a))') &
1781  ' DeltaIR1 = --------, DeltaIR2 = --------'
1782  write (50,'((a))') &
1783  ' epsIR epsIR^2'
1784  end if
1785  write (50,'((a))')
1786  write (50,'((a))') ' c(eps) = (4*pi)^eps\Gamma(1+eps), D = 4 -2*eps '
1787  write (50,'((a))')
1788  write (50,'((a))') ' you can freely choose the regularization parameters'
1789  if (rank.gt.0) then
1790  write (50,'((a))') ' of UV origin: muUV2 = mu^2, DeltaUV '
1791  end if
1792  write (50,'((a))') ' of IR origin: muIR2 = mu^2, DeltaIR1, DeltaIR2'
1793  write (50,'((a))')
1794  write (50,'((a))') ' note:'
1795  write (50,'((a))') ' - we effectively factor out a factor c(eps) '
1796  if (rank.gt.0) then
1797  write (50,'((a))') ' - by default DeltaUV = DeltaIR1 = DeltaIR2 = 0 '
1798  else
1799  write (50,'((a))') ' - by default DeltaIR1 = DeltaIR2 = 0 '
1800  end if
1801  write (50,'((a))') ' - suitable DeltaIR2 can be used to adapt the effective normalization'
1802 
1803  write (50,'((a))') ''
1804  write (50,'((a))') ' Results:'
1805 
1806  if (rank.gt.0) then
1807  do r=0,rank
1808  do n0=0,r/2
1809  do n1=0,r-n0
1810  do n2=0,r-n0-n1
1811  n3 = r-n0-n1-n2
1812  write (50,fmt10) n0,n1,n2,n3,cten(n0,n1,n2,n3)
1813  end do
1814  end do
1815  end do
1816  end do
1817 
1818  write (50,'(/(a))') ' Error estimates:'
1819  do r=0,rank
1820  write (50,fmt11) r,cerr(r)
1821  end do
1822  else
1823  write (50,fmt12) cten(0,0,0,0)
1824  end if
1825 
1826  write(*,'(/(a),(a)/)') ' The result has been written to the file ' &
1827  ,trim(fname)
1828 
1829  end subroutine writeresultcten
1830 
1831 
1832 
1833 
1834  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1835  ! subroutine writeresultD(caseN,casei,Dcoeff,Dcoeffuv,MomInv,masses2,rank,Derr)
1836  !
1837  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1838 
1839  subroutine writeresultd(caseN,casei,Dcoeff,Dcoeffuv,MomInv,masses2,rank,Derr)
1840  implicit none
1841 
1842  integer, intent(in) :: rank,caseN,casei
1843  double complex, intent(in) :: MomInv(6), masses2(0:3)
1844  double complex, intent(in) :: Dcoeff(0:rank/2,0:rank,0:rank,0:rank)
1845  double complex, intent(in) :: Dcoeffuv(0:rank/2,0:rank,0:rank,0:rank)
1846  double precision, optional, intent(in) :: Derr(0:rank)
1847  integer :: r,n0,n1,n2,n3,mode
1848  integer,parameter :: rankuv=4
1849  character(len=99) :: fname
1850  character(len=*),parameter :: fmt1 = "(A9,' = ',es23.16,' + i*',es23.16)"
1851  character(len=*),parameter :: fmt2 = "(A9,' = ',es23.16)"
1852  character(len=*),parameter :: fmt10 = "(' Dcoeff(',i1,3(',',i1),') = ',es23.16,' + i*',es23.16)"
1853  character(len=*),parameter :: fmt11 = "(' Derr(',i1,') = ',es23.16)"
1854  character(len=*),parameter :: fmt12 = "(' D0 = ',es23.16,' + i*',es23.16)"
1855 
1856 ! output file for result
1857  call getmode_cll(mode)
1858  select case (mode)
1859  case(1)
1860  fname = 'demo_4point_example00_coli.dat'
1861  case(2)
1862  fname = 'demo_4point_example00_dd.dat'
1863  case(3)
1864  fname = 'demo_4point_example00_comp.dat'
1865  end select
1866  write(fname(20:21),'(i2.2)') casei
1867 
1868 ! write(*,*) casei,fname
1869 
1870  open(unit=50,file=trim(fname),status='unknown')
1871 
1872  call getmuuv2_cll(muuv2)
1873  call getmuir2_cll(muir2)
1874  call getdeltauv_cll(deltauv)
1875  call getdeltair_cll(deltair1,deltair2)
1876 
1877  write (50,'(a37,i2,i3/)') ' Result for 4-point function, example',casen,casei
1878  write (50,'(a63,i2,i3,a)') ' The corresponding code can be found in demo.f90 under ''example',casen,casei,''''
1879  write (50,'((a))')
1880  write (50,'((a))') ' p31 '
1881  write (50,'((a))') ' ------------------ '
1882  write (50,'((a))') ' / \ '
1883  write (50,'((a))') ' m22 '
1884  write (50,'((a))') ' p21 --------------------- p32 \ '
1885  write (50,'((a))') ' | 2 | \ '
1886  write (50,'((a))') ' | | |'
1887  write (50,'((a))') ' m12 |1 3| m32 | p20'
1888  write (50,'((a))') ' | | |'
1889  write (50,'((a))') ' | 0 | / '
1890  write (50,'((a))') ' p10 --------------------- p30 / '
1891  write (50,'((a))') ' m02 '
1892  write (50,'((a))') ''
1893  write (50,'((a))') ' Input:'
1894  write (50,fmt1) ' p10 ',mominv(1)
1895  write (50,fmt1) ' p21 ',mominv(2)
1896  write (50,fmt1) ' p32 ',mominv(3)
1897  write (50,fmt1) ' p30 ',mominv(4)
1898  write (50,fmt1) ' p20 ',mominv(5)
1899  write (50,fmt1) ' p31 ',mominv(6)
1900  write (50,fmt1) ' m02 ',masses2(0)
1901  write (50,fmt1) ' m12 ',masses2(1)
1902  write (50,fmt1) ' m22 ',masses2(2)
1903  write (50,fmt1) ' m32 ',masses2(3)
1904  if (rank.ge.rankuv) then
1905  write (50,fmt2) ' muUV2 ',muuv2
1906  end if
1907  write (50,fmt2) ' muIR2 ',muir2
1908  if (rank.ge.rankuv) then
1909  write (50,fmt2) ' DeltaUV ',deltauv
1910  end if
1911  write (50,fmt2) ' DeltaIR1',deltair2
1912  write (50,fmt2) ' DeltaIR2',deltair1
1913  write (50,'((a))') ''
1914  write (50,'((a))') ' Conventions:'
1915  write (50,'((a))') ''
1916  if(rank.eq.0) then
1917  write (50,'((a))') ' (2*pi*mu)^(4-D) '
1918  write (50,'((a))') ' D0 = --------------- \int d^D q f(q,p_i)'
1919  write (50,'((a))') ' i*pi^2 '
1920 
1921  write (50,'((a))')
1922  write (50,'((a))') ' = D0_fin(muUV2,muIR2) + '// &
1923  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
1924  else if(rank.ge.rankuv) then
1925  write (50,'((a))') ' (2*pi*mu)^(4-D) '
1926  write (50,'((a))') ' D = --------------- \int d^D q f(q,p_i)'
1927  write (50,'((a))') ' i*pi^2 '
1928 
1929  write (50,'((a))')
1930  write (50,'((a))') ' = D_fin(muUV2,muIR2) + a_UV*DeltaUV + '// &
1931  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
1932  else
1933  write (50,'((a))') ' (2*pi*mu)^(4-D) '
1934  write (50,'((a))') ' D = --------------- \int d^D q f(q,p_i)'
1935  write (50,'((a))') ' i*pi^2 '
1936 
1937  write (50,'((a))')
1938  write (50,'((a))') ' = D_fin(muUV2,muIR2) + '// &
1939  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
1940  end if
1941  write (50,'((a))')
1942  write (50,'((a))') ' where'
1943  write (50,'((a))')
1944  if (rank.gt.rankuv) then
1945  write (50,'((a))') &
1946  ' c(epsUV) c(epsIR) c(epsIR)'
1947  write (50,'((a))') &
1948  ' DeltaUV = --------, DeltaIR1 = --------, DeltaIR2 = --------'
1949  write (50,'((a))') &
1950  ' epsUV epsIR epsIR^2'
1951  else
1952  write (50,'((a))') &
1953  ' c(epsIR) c(epsIR)'
1954  write (50,'((a))') &
1955  ' DeltaIR1 = --------, DeltaIR2 = --------'
1956  write (50,'((a))') &
1957  ' epsIR epsIR^2'
1958  end if
1959  write (50,'((a))')
1960  write (50,'((a))') ' c(eps) = (4*pi)^eps\Gamma(1+eps), D = 4 -2*eps '
1961  write (50,'((a))')
1962  write (50,'((a))') ' you can freely choose the regularization parameters'
1963  if (rank.gt.0) then
1964  write (50,'((a))') ' of UV origin: muUV2 = mu^2, DeltaUV '
1965  end if
1966  write (50,'((a))') ' of IR origin: muIR2 = mu^2, DeltaIR1, DeltaIR2'
1967  write (50,'((a))')
1968  write (50,'((a))') ' note:'
1969  write (50,'((a))') ' - we effectively factor out a factor c(eps) '
1970  if (rank.gt.0) then
1971  write (50,'((a))') ' - by default DeltaUV = DeltaIR1 = DeltaIR2 = 0 '
1972  else
1973  write (50,'((a))') ' - by default DeltaIR1 = DeltaIR2 = 0 '
1974  end if
1975  write (50,'((a))') ' - suitable DeltaIR2 can be used to adapt the effective normalization'
1976 
1977  write (50,'((a))') ''
1978  write (50,'((a))') ' Results:'
1979 
1980  if (rank.gt.0) then
1981  do r=0,rank
1982  do n0=0,r/2
1983  do n1=0,r-2*n0
1984  do n2=0,r-2*n0-n1
1985  n3 = r-2*n0-n1-n2
1986  write (50,fmt10) n0,n1,n2,n3,dcoeff(n0,n1,n2,n3)
1987  end do
1988  end do
1989  end do
1990  end do
1991 
1992  write (50,'(/(a))') ' Error estimates:'
1993  do r=0,rank
1994  write (50,fmt11) r,derr(r)
1995  end do
1996  else
1997  write (50,fmt12) dcoeff(0,0,0,0)
1998  end if
1999 
2000  write(*,'(/(a),(a)/)') ' The result has been written to the file ' &
2001  ,trim(fname)
2002 
2003  end subroutine writeresultd
2004 
2005 
2006  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2007  ! subroutine writeresultD0(caseN,casei,Dcoeff,Dcoeffuv,MomInv,masses2)
2008  !
2009  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2010 
2011  subroutine writeresultd0(caseN,casei,D0,MomInv,masses2)
2012  implicit none
2013 
2014  integer, intent(in) :: caseN,casei
2015  double complex, intent(in) :: D0
2016  double complex, intent(in) :: MomInv(6), masses2(0:3)
2017  double complex :: Dcoeff(0:0,0:0,0:0,0:0)
2018  double complex :: Dcoeffuv(0:0,0:0,0:0,0:0)
2019 
2020 
2021  dcoeff(0,0,0,0) = d0
2022  dcoeffuv(0,0,0,0) = 0d0
2023  call writeresultd(casen,casei,dcoeff,dcoeffuv,mominv,masses2,0)
2024 
2025  end subroutine writeresultd0
2026 
2027 
2028 
2029 
2030  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2031  ! subroutine writeresultDten(caseN,casei,Dten,Dtenuv,MomVec,MomInv,masses2,rank,Derr)
2032  !
2033  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2034 
2035  subroutine writeresultdten(caseN,casei,Dten,Dtenuv,MomVec,MomInv,masses2,rank,Derr)
2036  implicit none
2037 
2038  integer, intent(in) :: rank,caseN,casei
2039  double complex, intent(in) :: MomVec(0:3,3), MomInv(6), masses2(0:3)
2040  double complex, intent(in) :: Dten(0:rank,0:rank,0:rank,0:rank)
2041  double complex, intent(in) :: Dtenuv(0:rank,0:rank,0:rank,0:rank)
2042  double precision, optional, intent(in) :: Derr(0:rank)
2043  integer :: r,n0,n1,n2,n3,mode
2044  integer,parameter :: rankuv=4
2045  character(len=99) :: fname
2046  character(len=*),parameter :: fmt1 = "(A17,' = ',es23.16,' + i*',es23.16)"
2047  character(len=*),parameter :: fmt2 = "(A17,' = ',es23.16)"
2048  character(len=*),parameter :: fmt10 = "(' Dten(',i1,3(',',i1),') = ',es23.16,' + i*',es23.16)"
2049  character(len=*),parameter :: fmt11 = "(' Derr(',i1,') = ',es23.16)"
2050  character(len=*),parameter :: fmt12 = "(' D0 = ',es23.16,' + i*',es23.16)"
2051 
2052 ! output file for result
2053  call getmode_cll(mode)
2054  select case (mode)
2055  case(1)
2056  fname = 'demo_4point_example00_coli.dat'
2057  case(2)
2058  fname = 'demo_4point_example00_dd.dat'
2059  case(3)
2060  fname = 'demo_4point_example00_comp.dat'
2061  end select
2062  write(fname(20:21),'(i2.2)') casei
2063 
2064 ! write(*,*) casei,fname
2065 
2066  open(unit=50,file=trim(fname),status='unknown')
2067 
2068  call getmuuv2_cll(muuv2)
2069  call getmuir2_cll(muir2)
2070  call getdeltauv_cll(deltauv)
2071  call getdeltair_cll(deltair1,deltair2)
2072 
2073  write (50,'(a37,i2,i3/)') ' Result for 4-point function, example',casen,casei
2074  write (50,'(a63,i2,i3,a)') ' The corresponding code can be found in demo.f90 under ''example',casen,casei,''''
2075  write (50,'((a))')
2076  write (50,'((a))') ' p31 '
2077  write (50,'((a))') ' ------------------ '
2078  write (50,'((a))') ' / \ '
2079  write (50,'((a))') ' m22,p2vec '
2080  write (50,'((a))') ' p21 --------------------- p32 \ '
2081  write (50,'((a))') ' | 2 | \ '
2082  write (50,'((a))') ' | | |'
2083  write (50,'((a))') ' m12,p1vec |1 3| m32,p3vec | p20'
2084  write (50,'((a))') ' | | |'
2085  write (50,'((a))') ' | 0 | / '
2086  write (50,'((a))') ' p10 --------------------- p30 / '
2087  write (50,'((a))') ' m02 '
2088  write (50,'((a))') ''
2089  write (50,'((a))') ' Input:'
2090  write (50,fmt1) ' p1vec(0) ',momvec(0,1)
2091  write (50,fmt1) ' p1vec(1) ',momvec(1,1)
2092  write (50,fmt1) ' p1vec(2) ',momvec(2,1)
2093  write (50,fmt1) ' p1vec(3) ',momvec(3,1)
2094  write (50,fmt1) ' p2vec(0) ',momvec(0,2)
2095  write (50,fmt1) ' p2vec(1) ',momvec(1,2)
2096  write (50,fmt1) ' p2vec(2) ',momvec(2,2)
2097  write (50,fmt1) ' p2vec(3) ',momvec(3,2)
2098  write (50,fmt1) ' p3vec(0) ',momvec(0,3)
2099  write (50,fmt1) ' p3vec(1) ',momvec(1,3)
2100  write (50,fmt1) ' p3vec(2) ',momvec(2,3)
2101  write (50,fmt1) ' p3vec(3) ',momvec(3,3)
2102  write (50,fmt1) ' p10 ',mominv(1)
2103  write (50,fmt1) ' p21 ',mominv(2)
2104  write (50,fmt1) ' p32 ',mominv(3)
2105  write (50,fmt1) ' p30 ',mominv(4)
2106  write (50,fmt1) ' p20 ',mominv(5)
2107  write (50,fmt1) ' p31 ',mominv(6)
2108  write (50,fmt1) ' m02 ',masses2(0)
2109  write (50,fmt1) ' m12 ',masses2(1)
2110  write (50,fmt1) ' m22 ',masses2(2)
2111  write (50,fmt1) ' m32 ',masses2(3)
2112  if (rank.ge.rankuv) then
2113  write (50,fmt2) ' muUV2 ',muuv2
2114  end if
2115  write (50,fmt2) ' muIR2 ',muir2
2116  if (rank.ge.rankuv) then
2117  write (50,fmt2) ' DeltaUV ',deltauv
2118  end if
2119  write (50,fmt2) ' DeltaIR1 ',deltair2
2120  write (50,fmt2) ' DeltaIR2 ',deltair1
2121  write (50,'((a))') ''
2122  write (50,'((a))') ' Conventions:'
2123  write (50,'((a))') ''
2124  if(rank.eq.0) then
2125  write (50,'((a))') ' (2*pi*mu)^(4-D) '
2126  write (50,'((a))') ' D0 = --------------- \int d^D q f(q,p_i)'
2127  write (50,'((a))') ' i*pi^2 '
2128 
2129  write (50,'((a))')
2130  write (50,'((a))') ' = D0_fin(muUV2,muIR2) + '// &
2131  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
2132  else if(rank.ge.rankuv) then
2133  write (50,'((a))') ' (2*pi*mu)^(4-D) '
2134  write (50,'((a))') ' D = --------------- \int d^D q f(q,p_i)'
2135  write (50,'((a))') ' i*pi^2 '
2136 
2137  write (50,'((a))')
2138  write (50,'((a))') ' = D_fin(muUV2,muIR2) + a_UV*DeltaUV + '// &
2139  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
2140  else
2141  write (50,'((a))') ' (2*pi*mu)^(4-D) '
2142  write (50,'((a))') ' D = --------------- \int d^D q f(q,p_i)'
2143  write (50,'((a))') ' i*pi^2 '
2144 
2145  write (50,'((a))')
2146  write (50,'((a))') ' = D_fin(muUV2,muIR2) + '// &
2147  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
2148  end if
2149  write (50,'((a))')
2150  write (50,'((a))') ' where'
2151  write (50,'((a))')
2152  if (rank.gt.rankuv) then
2153  write (50,'((a))') &
2154  ' c(epsUV) c(epsIR) c(epsIR)'
2155  write (50,'((a))') &
2156  ' DeltaUV = --------, DeltaIR1 = --------, DeltaIR2 = --------'
2157  write (50,'((a))') &
2158  ' epsUV epsIR epsIR^2'
2159  else
2160  write (50,'((a))') &
2161  ' c(epsIR) c(epsIR)'
2162  write (50,'((a))') &
2163  ' DeltaIR1 = --------, DeltaIR2 = --------'
2164  write (50,'((a))') &
2165  ' epsIR epsIR^2'
2166  end if
2167  write (50,'((a))')
2168  write (50,'((a))') ' c(eps) = (4*pi)^eps\Gamma(1+eps), D = 4 -2*eps '
2169  write (50,'((a))')
2170  write (50,'((a))') ' you can freely choose the regularization parameters'
2171  if (rank.gt.0) then
2172  write (50,'((a))') ' of UV origin: muUV2 = mu^2, DeltaUV '
2173  end if
2174  write (50,'((a))') ' of IR origin: muIR2 = mu^2, DeltaIR1, DeltaIR2'
2175  write (50,'((a))')
2176  write (50,'((a))') ' note:'
2177  write (50,'((a))') ' - we effectively factor out a factor c(eps) '
2178  if (rank.gt.0) then
2179  write (50,'((a))') ' - by default DeltaUV = DeltaIR1 = DeltaIR2 = 0 '
2180  else
2181  write (50,'((a))') ' - by default DeltaIR1 = DeltaIR2 = 0 '
2182  end if
2183  write (50,'((a))') ' - suitable DeltaIR2 can be used to adapt the effective normalization'
2184 
2185  write (50,'((a))') ''
2186  write (50,'((a))') ' Results:'
2187 
2188  if (rank.gt.0) then
2189  do r=0,rank
2190  do n0=0,r
2191  do n1=0,r-n0
2192  do n2=0,r-n0-n1
2193  n3 = r-n0-n1-n2
2194  write (50,fmt10) n0,n1,n2,n3,dten(n0,n1,n2,n3)
2195  end do
2196  end do
2197  end do
2198  end do
2199 
2200  write (50,'(/(a))') ' Error estimates:'
2201  do r=0,rank
2202  write (50,fmt11) r,derr(r)
2203  end do
2204  else
2205  write (50,fmt12) dten(0,0,0,0)
2206  end if
2207 
2208  write(*,'(/(a),(a)/)') ' The result has been written to the file ' &
2209  ,trim(fname)
2210 
2211  end subroutine writeresultdten
2212 
2213 
2214 
2215 
2216  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2217  ! subroutine writeresultE(caseN,casei,Ecoeff,Ecoeffuv,MomInv,masses2,rank,Eerr)
2218  !
2219  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2220 
2221  subroutine writeresulte(caseN,casei,Ecoeff,Ecoeffuv,MomInv,masses2,rank,Eerr)
2222  implicit none
2223 
2224  integer, intent(in) :: rank,caseN,casei
2225  double complex, intent(in) :: MomInv(10), masses2(0:4)
2226  double complex, intent(in) :: Ecoeff(0:rank/2,0:rank,0:rank,0:rank,0:rank)
2227  double complex, intent(in) :: Ecoeffuv(0:rank/2,0:rank,0:rank,0:rank,0:rank)
2228  double precision, optional, intent(in) :: Eerr(0:rank)
2229  integer :: r,n0,n1,n2,n3,n4,mode
2230  integer,parameter :: rankuv=6
2231  character(len=99) :: fname
2232  character(len=*),parameter :: fmt1 = "(A9,' = ',es23.16,' + i*',es23.16)"
2233  character(len=*),parameter :: fmt2 = "(A9,' = ',es23.16)"
2234  character(len=*),parameter :: fmt10 = "(' Ecoeff(',i1,4(',',i1),') = ',es23.16,' + i*',es23.16)"
2235  character(len=*),parameter :: fmt11 = "(' Eerr(',i1,') = ',es23.16)"
2236  character(len=*),parameter :: fmt12 = "(' E0 = ',es23.16,' + i*',es23.16)"
2237 
2238 ! output file for result
2239  call getmode_cll(mode)
2240  select case (mode)
2241  case(1)
2242  fname = 'demo_5point_example00_coli.dat'
2243  case(2)
2244  fname = 'demo_5point_example00_dd.dat'
2245  case(3)
2246  fname = 'demo_5point_example00_comp.dat'
2247  end select
2248  write(fname(20:21),'(i2.2)') casei
2249 
2250 ! write(*,*) casei,fname
2251 
2252  open(unit=50,file=trim(fname),status='unknown')
2253 
2254  call getmuuv2_cll(muuv2)
2255  call getmuir2_cll(muir2)
2256  call getdeltauv_cll(deltauv)
2257  call getdeltair_cll(deltair1,deltair2)
2258 
2259  write (50,'(a37,i2,i3/)') ' Result for 5-point function, example',casen,casei
2260  write (50,'(a63,i2,i3,a)') ' The corresponding code can be found in demo.f90 under ''example',casen,casei,''''
2261  write (50,'((a))')
2262  write (50,'((a))') ' p31 '
2263  write (50,'((a))') ' ---------^--------- '
2264  write (50,'((a))') ' / \ '
2265  write (50,'((a))') ' m22 '
2266  write (50,'((a))') ' / p21 --------------------- p32 \ '
2267  write (50,'((a))') ' / | 2 \ \ '
2268  write (50,'((a))') ' / | 3\ m32 | p42'
2269  write (50,'((a))') ' | |1 \ / '
2270  write (50,'((a))') ' p20 | m12 | >---- p43 < '
2271  write (50,'((a))') ' | | / \ '
2272  write (50,'((a))') ' \ | 4/ m42 | p30'
2273  write (50,'((a))') ' \ | 0 / / '
2274  write (50,'((a))') ' \ p10 --------------------- p40 / '
2275  write (50,'((a))') ' m02 '
2276  write (50,'((a))') ' \ / '
2277  write (50,'((a))') ' ---------v--------- '
2278  write (50,'((a))') ' p41 '
2279  write (50,'((a))') ''
2280  write (50,'((a))') ' Input:'
2281  write (50,fmt1) ' p10 ',mominv(1)
2282  write (50,fmt1) ' p21 ',mominv(2)
2283  write (50,fmt1) ' p32 ',mominv(3)
2284  write (50,fmt1) ' p43 ',mominv(4)
2285  write (50,fmt1) ' p40 ',mominv(5)
2286  write (50,fmt1) ' p20 ',mominv(6)
2287  write (50,fmt1) ' p31 ',mominv(7)
2288  write (50,fmt1) ' p42 ',mominv(8)
2289  write (50,fmt1) ' p30 ',mominv(9)
2290  write (50,fmt1) ' p41 ',mominv(10)
2291  write (50,fmt1) ' m02 ',masses2(0)
2292  write (50,fmt1) ' m12 ',masses2(1)
2293  write (50,fmt1) ' m22 ',masses2(2)
2294  write (50,fmt1) ' m32 ',masses2(3)
2295  write (50,fmt1) ' m42 ',masses2(4)
2296  if (rank.ge.rankuv) then
2297  write (50,fmt2) ' muUV2 ',muuv2
2298  end if
2299  write (50,fmt2) ' muIR2 ',muir2
2300  if (rank.ge.rankuv) then
2301  write (50,fmt2) ' deltaUV ',deltauv
2302  end if
2303  write (50,fmt2) ' DeltaIR1',deltair2
2304  write (50,fmt2) ' DeltaIR2',deltair1
2305  write (50,'((a))') ''
2306  write (50,'((a))') ' Conventions:'
2307  write (50,'((a))') ''
2308  if(rank.eq.0) then
2309  write (50,'((a))') ' (2*pi*mu)^(4-D) '
2310  write (50,'((a))') ' E0 = --------------- \int d^D q f(q,p_i)'
2311  write (50,'((a))') ' i*pi^2 '
2312 
2313  write (50,'((a))')
2314  write (50,'((a))') ' = E0_fin(muUV2,muIR2) + '// &
2315  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
2316  elseif(rank.ge.rankuv) then
2317  write (50,'((a))') ' (2*pi*mu)^(4-D) '
2318  write (50,'((a))') ' E = --------------- \int d^D q f(q,p_i)'
2319  write (50,'((a))') ' i*pi^2 '
2320 
2321  write (50,'((a))')
2322  write (50,'((a))') ' = E_fin(muUV2,muIR2) + a_UV*DeltaUV + '// &
2323  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
2324  else
2325  write (50,'((a))') ' (2*pi*mu)^(4-D) '
2326  write (50,'((a))') ' E = --------------- \int d^D q f(q,p_i)'
2327  write (50,'((a))') ' i*pi^2 '
2328 
2329  write (50,'((a))')
2330  write (50,'((a))') ' = E_fin(muUV2,muIR2) + '// &
2331  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
2332  end if
2333  write (50,'((a))')
2334  write (50,'((a))') ' where'
2335  write (50,'((a))')
2336  if (rank.ge.rankuv) then
2337  write (50,'((a))') &
2338  ' c(epsUV) c(epsIR) c(epsIR)'
2339  write (50,'((a))') &
2340  ' DeltaUV = --------, DeltaIR1 = --------, DeltaIR2 = --------'
2341  write (50,'((a))') &
2342  ' epsUV epsIR epsIR^2'
2343  else
2344  write (50,'((a))') &
2345  ' c(epsIR) c(epsIR)'
2346  write (50,'((a))') &
2347  ' DeltaIR1 = --------, DeltaIR2 = --------'
2348  write (50,'((a))') &
2349  ' epsIR epsIR^2'
2350  end if
2351  write (50,'((a))')
2352  write (50,'((a))') ' c(eps) = (4*pi)^eps\Gamma(1+eps), D = 4 -2*eps '
2353  write (50,'((a))')
2354  write (50,'((a))') ' you can freely choose the regularization parameters'
2355  if (rank.ge.rankuv) then
2356  write (50,'((a))') ' of UV origin: muUV2 = mu^2, DeltaUV '
2357  end if
2358  write (50,'((a))') ' of IR origin: muIR2 = mu^2, DeltaIR1, DeltaIR2'
2359  write (50,'((a))')
2360  write (50,'((a))') ' note:'
2361  write (50,'((a))') ' - we effectively factor out a factor c(eps) '
2362  if (rank.ge.rankuv) then
2363  write (50,'((a))') ' - by default DeltaUV = DeltaIR1 = DeltaIR2 = 0 '
2364  else
2365  write (50,'((a))') ' - by default DeltaIR1 = DeltaIR2 = 0 '
2366  end if
2367  write (50,'((a))') ' - suitable DeltaIR2 can be used to adapt the effective normalization'
2368 
2369  write (50,'((a))') ''
2370  write (50,'((a))') ' Results:'
2371 
2372  if (rank.gt.0) then
2373  do r=0,rank
2374  do n0=0,r/2
2375  do n1=0,r-2*n0
2376  do n2=0,r-2*n0-n1
2377  do n3=0,r-2*n0-n1-n2
2378  n4 = r-2*n0-n1-n2-n3
2379  write (50,fmt10) n0,n1,n2,n3,n4,ecoeff(n0,n1,n2,n3,n4)
2380  end do
2381  end do
2382  end do
2383  end do
2384  end do
2385 
2386  write (50,'(/(a))') ' Error estimates:'
2387  do r=0,rank
2388  write (50,fmt11) r,eerr(r)
2389  end do
2390  else
2391  write (50,fmt12) ecoeff(0,0,0,0,0)
2392  end if
2393 
2394  write(*,'(/(a),(a)/)') ' The result has been written to the file ' &
2395  ,trim(fname)
2396 
2397  end subroutine writeresulte
2398 
2399 
2400  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2401  ! subroutine writeresultE0(caseN,casei,Ecoeff,Ecoeffuv,MomInv,masses2)
2402  !
2403  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2404 
2405  subroutine writeresulte0(caseN,casei,E0,MomInv,masses2)
2406  implicit none
2407 
2408  integer, intent(in) :: caseN,casei
2409  double complex, intent(in) :: E0
2410  double complex, intent(in) :: MomInv(10), masses2(0:4)
2411  double complex :: Ecoeff(0:0,0:0,0:0,0:0,0:0)
2412  double complex :: Ecoeffuv(0:0,0:0,0:0,0:0,0:0)
2413 
2414 
2415  ecoeff(0,0,0,0,0) = e0
2416  ecoeffuv(0,0,0,0,0) = 0d0
2417  call writeresulte(casen,casei,ecoeff,ecoeffuv,mominv,masses2,0)
2418 
2419  end subroutine writeresulte0
2420 
2421 
2422 
2423 
2424  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2425  ! subroutine writeresultEten(caseN,casei,Eten,Etenuv,MomVec,MomInv,masses2,rank,Eerr)
2426  !
2427  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2428 
2429  subroutine writeresulteten(caseN,casei,Eten,Etenuv,MomVec,MomInv,masses2,rank,Eerr)
2430  implicit none
2431 
2432  integer, intent(in) :: rank,caseN,casei
2433  double complex, intent(in) :: MomVec(0:3,4), MomInv(10), masses2(0:4)
2434  double complex, intent(in) :: Eten(0:rank,0:rank,0:rank,0:rank)
2435  double complex, intent(in) :: Etenuv(0:rank,0:rank,0:rank,0:rank)
2436  double precision, optional, intent(in) :: Eerr(0:rank)
2437  integer :: r,n0,n1,n2,n3,mode
2438  integer,parameter :: rankuv=6
2439  character(len=99) :: fname
2440  character(len=*),parameter :: fmt1 = "(A17,' = ',es23.16,' + i*',es23.16)"
2441  character(len=*),parameter :: fmt2 = "(A17,' = ',es23.16)"
2442  character(len=*),parameter :: fmt10 = "(' Eten(',i1,3(',',i1),') = ',es23.16,' + i*',es23.16)"
2443  character(len=*),parameter :: fmt11 = "(' Eerr(',i1,') = ',es23.16)"
2444  character(len=*),parameter :: fmt12 = "(' E0 = ',es23.16,' + i*',es23.16)"
2445 
2446 ! output file for result
2447  call getmode_cll(mode)
2448  select case (mode)
2449  case(1)
2450  fname = 'demo_5point_example00_coli.dat'
2451  case(2)
2452  fname = 'demo_5point_example00_dd.dat'
2453  case(3)
2454  fname = 'demo_5point_example00_comp.dat'
2455  end select
2456  write(fname(20:21),'(i2.2)') casei
2457 
2458 ! write(*,*) casei,fname
2459 
2460  open(unit=50,file=trim(fname),status='unknown')
2461 
2462  call getmuuv2_cll(muuv2)
2463  call getmuir2_cll(muir2)
2464  call getdeltauv_cll(deltauv)
2465  call getdeltair_cll(deltair1,deltair2)
2466 
2467  write (50,'(a37,i2,i3/)') ' Result for 5-point function, example',casen,casei
2468  write (50,'(a63,i2,i3,a)') ' The corresponding code can be found in demo.f90 under ''example',casen,casei,''''
2469  write (50,'((a))')
2470  write (50,'((a))') ' p31 '
2471  write (50,'((a))') ' ---------^--------- '
2472  write (50,'((a))') ' / \ '
2473  write (50,'((a))') ' m22,p2vec '
2474  write (50,'((a))') ' / p21 --------------------- p32 \ '
2475  write (50,'((a))') ' / | 2 \ \ '
2476  write (50,'((a))') ' / | 3\ m32,p3vec | p42'
2477  write (50,'((a))') ' | |1 \ / '
2478  write (50,'((a))') ' p20 | m12,p1vec| >---- p43 < '
2479  write (50,'((a))') ' | | / \ '
2480  write (50,'((a))') ' \ | 4/ m42,p4vec | p30'
2481  write (50,'((a))') ' \ | 0 / / '
2482  write (50,'((a))') ' \ p10 --------------------- p40 / '
2483  write (50,'((a))') ' m02 '
2484  write (50,'((a))') ' \ / '
2485  write (50,'((a))') ' ---------v--------- '
2486  write (50,'((a))') ' p41 '
2487  write (50,'((a))') ''
2488  write (50,'((a))') ' Input:'
2489  write (50,fmt1) ' p1vec(0) ',momvec(0,1)
2490  write (50,fmt1) ' p1vec(1) ',momvec(1,1)
2491  write (50,fmt1) ' p1vec(2) ',momvec(2,1)
2492  write (50,fmt1) ' p1vec(3) ',momvec(3,1)
2493  write (50,fmt1) ' p2vec(0) ',momvec(0,2)
2494  write (50,fmt1) ' p2vec(1) ',momvec(1,2)
2495  write (50,fmt1) ' p2vec(2) ',momvec(2,2)
2496  write (50,fmt1) ' p2vec(3) ',momvec(3,2)
2497  write (50,fmt1) ' p3vec(0) ',momvec(0,3)
2498  write (50,fmt1) ' p3vec(1) ',momvec(1,3)
2499  write (50,fmt1) ' p3vec(2) ',momvec(2,3)
2500  write (50,fmt1) ' p3vec(3) ',momvec(3,3)
2501  write (50,fmt1) ' p4vec(0) ',momvec(0,4)
2502  write (50,fmt1) ' p4vec(1) ',momvec(1,4)
2503  write (50,fmt1) ' p4vec(2) ',momvec(2,4)
2504  write (50,fmt1) ' p4vec(3) ',momvec(3,4)
2505  write (50,fmt1) ' p10 ',mominv(1)
2506  write (50,fmt1) ' p21 ',mominv(2)
2507  write (50,fmt1) ' p32 ',mominv(3)
2508  write (50,fmt1) ' p43 ',mominv(4)
2509  write (50,fmt1) ' p40 ',mominv(5)
2510  write (50,fmt1) ' p20 ',mominv(6)
2511  write (50,fmt1) ' p31 ',mominv(7)
2512  write (50,fmt1) ' p42 ',mominv(8)
2513  write (50,fmt1) ' p30 ',mominv(9)
2514  write (50,fmt1) ' p41 ',mominv(10)
2515  write (50,fmt1) ' m02 ',masses2(0)
2516  write (50,fmt1) ' m12 ',masses2(1)
2517  write (50,fmt1) ' m22 ',masses2(2)
2518  write (50,fmt1) ' m32 ',masses2(3)
2519  write (50,fmt1) ' m42 ',masses2(4)
2520  if (rank.ge.rankuv) then
2521  write (50,fmt2) ' muUV2 ',muuv2
2522  end if
2523  write (50,fmt2) ' muIR2 ',muir2
2524  if (rank.ge.rankuv) then
2525  write (50,fmt2) ' deltaUV ',deltauv
2526  end if
2527  write (50,fmt2) ' DeltaIR1 ',deltair2
2528  write (50,fmt2) ' DeltaIR2 ',deltair1
2529  write (50,'((a))') ''
2530  write (50,'((a))') ' Conventions:'
2531  write (50,'((a))') ''
2532  if(rank.eq.0) then
2533  write (50,'((a))') ' (2*pi*mu)^(4-D) '
2534  write (50,'((a))') ' E0 = --------------- \int d^D q f(q,p_i)'
2535  write (50,'((a))') ' i*pi^2 '
2536 
2537  write (50,'((a))')
2538  write (50,'((a))') ' = E0_fin(muUV2,muIR2) + '// &
2539  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
2540  elseif(rank.ge.rankuv) then
2541  write (50,'((a))') ' (2*pi*mu)^(4-D) '
2542  write (50,'((a))') ' E = --------------- \int d^D q f(q,p_i)'
2543  write (50,'((a))') ' i*pi^2 '
2544 
2545  write (50,'((a))')
2546  write (50,'((a))') ' = E_fin(muUV2,muIR2) + a_UV*DeltaUV + '// &
2547  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
2548  else
2549  write (50,'((a))') ' (2*pi*mu)^(4-D) '
2550  write (50,'((a))') ' E = --------------- \int d^D q f(q,p_i)'
2551  write (50,'((a))') ' i*pi^2 '
2552 
2553  write (50,'((a))')
2554  write (50,'((a))') ' = E_fin(muUV2,muIR2) + '// &
2555  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
2556  end if
2557  write (50,'((a))')
2558  write (50,'((a))') ' where'
2559  write (50,'((a))')
2560  if (rank.ge.rankuv) then
2561  write (50,'((a))') &
2562  ' c(epsUV) c(epsIR) c(epsIR)'
2563  write (50,'((a))') &
2564  ' DeltaUV = --------, DeltaIR1 = --------, DeltaIR2 = --------'
2565  write (50,'((a))') &
2566  ' epsUV epsIR epsIR^2'
2567  else
2568  write (50,'((a))') &
2569  ' c(epsIR) c(epsIR)'
2570  write (50,'((a))') &
2571  ' DeltaIR1 = --------, DeltaIR2 = --------'
2572  write (50,'((a))') &
2573  ' epsIR epsIR^2'
2574  end if
2575  write (50,'((a))')
2576  write (50,'((a))') ' c(eps) = (4*pi)^eps\Gamma(1+eps), D = 4 -2*eps '
2577  write (50,'((a))')
2578  write (50,'((a))') ' you can freely choose the regularization parameters'
2579  if (rank.ge.rankuv) then
2580  write (50,'((a))') ' of UV origin: muUV2 = mu^2, DeltaUV '
2581  end if
2582  write (50,'((a))') ' of IR origin: muIR2 = mu^2, DeltaIR1, DeltaIR2'
2583  write (50,'((a))')
2584  write (50,'((a))') ' note:'
2585  write (50,'((a))') ' - we effectively factor out a factor c(eps) '
2586  if (rank.ge.rankuv) then
2587  write (50,'((a))') ' - by default DeltaUV = DeltaIR1 = DeltaIR2 = 0 '
2588  else
2589  write (50,'((a))') ' - by default DeltaIR1 = DeltaIR2 = 0 '
2590  end if
2591  write (50,'((a))') ' - suitable DeltaIR2 can be used to adapt the effective normalization'
2592 
2593  write (50,'((a))') ''
2594  write (50,'((a))') ' Results:'
2595 
2596  if (rank.gt.0) then
2597  do r=0,rank
2598  do n0=0,r
2599  do n1=0,r-n0
2600  do n2=0,r-n0-n1
2601  n3=r-n0-n1-n2
2602  write (50,fmt10) n0,n1,n2,n3,eten(n0,n1,n2,n3)
2603  end do
2604  end do
2605  end do
2606  end do
2607 
2608  write (50,'(/(a))') ' Error estimates:'
2609  do r=0,rank
2610  write (50,fmt11) r,eerr(r)
2611  end do
2612  else
2613  write (50,fmt12) eten(0,0,0,0)
2614  end if
2615 
2616  write(*,'(/(a),(a)/)') ' The result has been written to the file ' &
2617  ,trim(fname)
2618 
2619  end subroutine writeresulteten
2620 
2621 
2622 
2623 
2624 
2625 
2626  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2627  ! subroutine writeresultF(caseN,casei,Fcoeff,Fcoeffuv,MomInv,masses2,rank,Ferr)
2628  !
2629  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2630 
2631  subroutine writeresultf(caseN,casei,Fcoeff,Fcoeffuv,MomInv,masses2,rank,Ferr)
2632  implicit none
2633 
2634  integer, intent(in) :: rank,caseN,casei
2635  double complex, intent(in) :: MomInv(15), masses2(0:5)
2636  double complex, intent(in) :: Fcoeff(0:rank/2,0:rank,0:rank,0:rank,0:rank,0:rank)
2637  double complex, intent(in) :: Fcoeffuv(0:rank/2,0:rank,0:rank,0:rank,0:rank,0:rank)
2638  double precision, optional, intent(in) :: Ferr(0:rank)
2639  integer :: r,n0,n1,n2,n3,n4,n5,mode
2640  integer,parameter :: rankuv=8
2641  character(len=99) :: fname
2642  character(len=*),parameter :: fmt1 = "(A9,' = ',es23.16,' + i*',es23.16)"
2643  character(len=*),parameter :: fmt2 = "(A9,' = ',es23.16)"
2644  character(len=*),parameter :: fmt10 = "(' Fcoeff(',i1,5(',',i1),') = ',es23.16,' + i*',es23.16)"
2645  character(len=*),parameter :: fmt11 = "(' Ferr(',i1,') = ',es23.16)"
2646  character(len=*),parameter :: fmt12 = "(' F0 = ',es23.16,' + i*',es23.16)"
2647 
2648 ! output file for result
2649  call getmode_cll(mode)
2650  select case (mode)
2651  case(1)
2652  fname = 'demo_6point_example00_coli.dat'
2653  case(2)
2654  fname = 'demo_6point_example00_dd.dat'
2655  case(3)
2656  fname = 'demo_6point_example00_comp.dat'
2657  end select
2658  write(fname(20:21),'(i2.2)') casei
2659 
2660 ! write(*,*) casei,fname
2661 
2662  open(unit=50,file=trim(fname),status='unknown')
2663 
2664  call getmuuv2_cll(muuv2)
2665  call getmuir2_cll(muir2)
2666  call getdeltauv_cll(deltauv)
2667  call getdeltair_cll(deltair1,deltair2)
2668 
2669  write (50,'(a37,i2,i3/)') ' Result for 6-point function, example',casen,casei
2670  write (50,'(a63,i2,i3,a)') ' The corresponding code can be found in demo.f90 under ''example',casen,casei,''''
2671  write (50,'((a))')
2672  write (50,'((a))') ' p42 '
2673  write (50,'((a))') ' ---------^--------- '
2674  write (50,'((a))') ' / \ '
2675  write (50,'((a))') ' m32 '
2676  write (50,'((a))') ' / p32 --------------------- p43 \ '
2677  write (50,'((a))') ' / / 3 \ \ '
2678  write (50,'((a))') ' p31 | m22 /2 4\ m42 | p53'
2679  write (50,'((a))') ' \ / \ / '
2680  write (50,'((a))') ' > p21 ----< >---- p54 < '
2681  write (50,'((a))') ' / \ / \ '
2682  write (50,'((a))') ' p20 | m12 \1 5/ m52 | p40'
2683  write (50,'((a))') ' \ \ 0 / / '
2684  write (50,'((a))') ' \ p10 --------------------- p50 / '
2685  write (50,'((a))') ' m02 '
2686  write (50,'((a))') ' \ / '
2687  write (50,'((a))') ' ---------v--------- '
2688  write (50,'((a))') ' p51 '
2689  write (50,'((a))') ''
2690  write (50,'((a))') ' Input:'
2691  write (50,fmt1) ' p10 ',mominv(1)
2692  write (50,fmt1) ' p21 ',mominv(2)
2693  write (50,fmt1) ' p32 ',mominv(3)
2694  write (50,fmt1) ' p43 ',mominv(4)
2695  write (50,fmt1) ' p54 ',mominv(5)
2696  write (50,fmt1) ' p50 ',mominv(6)
2697  write (50,fmt1) ' p20 ',mominv(7)
2698  write (50,fmt1) ' p31 ',mominv(8)
2699  write (50,fmt1) ' p42 ',mominv(9)
2700  write (50,fmt1) ' p53 ',mominv(10)
2701  write (50,fmt1) ' p40 ',mominv(11)
2702  write (50,fmt1) ' p51 ',mominv(12)
2703  write (50,fmt1) ' p30 ',mominv(13)
2704  write (50,fmt1) ' p41 ',mominv(14)
2705  write (50,fmt1) ' p52 ',mominv(15)
2706  write (50,fmt1) ' m02 ',masses2(0)
2707  write (50,fmt1) ' m12 ',masses2(1)
2708  write (50,fmt1) ' m22 ',masses2(2)
2709  write (50,fmt1) ' m32 ',masses2(3)
2710  write (50,fmt1) ' m42 ',masses2(4)
2711  write (50,fmt1) ' m52 ',masses2(5)
2712  if (rank.ge.rankuv) then
2713  write (50,fmt2) ' muUV2 ',muuv2
2714  end if
2715  write (50,fmt2) ' muIR2 ',muir2
2716  if (rank.ge.rankuv) then
2717  write (50,fmt2) ' deltaUV ',deltauv
2718  end if
2719  write (50,fmt2) ' DeltaIR1',deltair2
2720  write (50,fmt2) ' DeltaIR2',deltair1
2721  write (50,'((a))') ''
2722  write (50,'((a))') ' Conventions:'
2723  write (50,'((a))') ''
2724  if(rank.eq.0) then
2725  write (50,'((a))') ' (2*pi*mu)^(4-D) '
2726  write (50,'((a))') ' F0 = --------------- \int d^D q f(q,p_i)'
2727  write (50,'((a))') ' i*pi^2 '
2728 
2729  write (50,'((a))')
2730  write (50,'((a))') ' = F0_fin(muUV2,muIR2) + '// &
2731  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
2732  elseif(rank.ge.rankuv) then
2733  write (50,'((a))') ' (2*pi*mu)^(4-D) '
2734  write (50,'((a))') ' F = --------------- \int d^D q f(q,p_i)'
2735  write (50,'((a))') ' i*pi^2 '
2736 
2737  write (50,'((a))')
2738  write (50,'((a))') ' = F_fin(muUV2,muIR2) + a_UV*DeltaUV + '// &
2739  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
2740  else
2741  write (50,'((a))') ' (2*pi*mu)^(4-D) '
2742  write (50,'((a))') ' F = --------------- \int d^D q f(q,p_i)'
2743  write (50,'((a))') ' i*pi^2 '
2744 
2745  write (50,'((a))')
2746  write (50,'((a))') ' = F_fin(muUV2,muIR2) + '// &
2747  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
2748  end if
2749  write (50,'((a))')
2750  write (50,'((a))') ' where'
2751  write (50,'((a))')
2752  if (rank.ge.rankuv) then
2753  write (50,'((a))') &
2754  ' c(epsUV) c(epsIR) c(epsIR)'
2755  write (50,'((a))') &
2756  ' DeltaUV = --------, DeltaIR1 = --------, DeltaIR2 = --------'
2757  write (50,'((a))') &
2758  ' epsUV epsIR epsIR^2'
2759  else
2760  write (50,'((a))') &
2761  ' c(epsIR) c(epsIR)'
2762  write (50,'((a))') &
2763  ' DeltaIR1 = --------, DeltaIR2 = --------'
2764  write (50,'((a))') &
2765  ' epsIR epsIR^2'
2766  end if
2767  write (50,'((a))')
2768  write (50,'((a))') ' c(eps) = (4*pi)^eps\Gamma(1+eps), D = 4 -2*eps '
2769  write (50,'((a))')
2770  write (50,'((a))') ' you can freely choose the regularization parameters'
2771  if (rank.ge.rankuv) then
2772  write (50,'((a))') ' of UV origin: muUV2 = mu^2, DeltaUV '
2773  end if
2774  write (50,'((a))') ' of IR origin: muIR2 = mu^2, DeltaIR1, DeltaIR2'
2775  write (50,'((a))')
2776  write (50,'((a))') ' note:'
2777  write (50,'((a))') ' - we effectively factor out a factor c(eps) '
2778  if (rank.ge.rankuv) then
2779  write (50,'((a))') ' - by default DeltaUV = DeltaIR1 = DeltaIR2 = 0 '
2780  else
2781  write (50,'((a))') ' - by default DeltaIR1 = DeltaIR2 = 0 '
2782  end if
2783  write (50,'((a))') ' - suitable DeltaIR2 can be used to adapt the effective normalization'
2784 
2785  write (50,'((a))') ''
2786  write (50,'((a))') ' Results:'
2787 
2788  if (rank.gt.0) then
2789  do r=0,rank
2790  do n0=0,r/2
2791  do n1=0,r-2*n0
2792  do n2=0,r-2*n0-n1
2793  do n3=0,r-2*n0-n1-n2
2794  do n4=0,r-2*n0-n1-n2-n3
2795  n5 = r-2*n0-n1-n2-n3-n4
2796  write (50,fmt10) n0,n1,n2,n3,n4,n5,fcoeff(n0,n1,n2,n3,n4,n5)
2797  end do
2798  end do
2799  end do
2800  end do
2801  end do
2802  end do
2803 
2804  write (50,'(/(a))') ' Error estimates:'
2805  do r=0,rank
2806  write (50,fmt11) r,ferr(r)
2807  end do
2808  else
2809  write (50,fmt12) fcoeff(0,0,0,0,0,0)
2810  end if
2811 
2812  write(*,'(/(a),(a)/)') ' The result has been written to the file ' &
2813  ,trim(fname)
2814 
2815  end subroutine writeresultf
2816 
2817 
2818  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2819  ! subroutine writeresultF0(caseN,casei,Fcoeff,Fcoeffuv,MomInv,masses2)
2820  !
2821  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2822 
2823  subroutine writeresultf0(caseN,casei,F0,MomInv,masses2)
2824  implicit none
2825 
2826  integer, intent(in) :: caseN,casei
2827  double complex, intent(in) :: F0
2828  double complex, intent(in) :: MomInv(15), masses2(0:5)
2829  double complex :: Fcoeff(0:0,0:0,0:0,0:0,0:0,0:0)
2830  double complex :: Fcoeffuv(0:0,0:0,0:0,0:0,0:0,0:0)
2831 
2832 
2833  fcoeff(0,0,0,0,0,0) = f0
2834  fcoeffuv(0,0,0,0,0,0) = 0d0
2835  call writeresultf(casen,casei,fcoeff,fcoeffuv,mominv,masses2,0)
2836 
2837  end subroutine writeresultf0
2838 
2839 
2840 
2841 
2842 
2843 
2844  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2845  ! subroutine writeresultFten(caseN,casei,Ften,Ftenuv,MomVec,MomInv,masses2,rank,Ferr)
2846  !
2847  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2848 
2849  subroutine writeresultften(caseN,casei,Ften,Ftenuv,MomVec,MomInv,masses2,rank,Ferr)
2850  implicit none
2851 
2852  integer, intent(in) :: rank,caseN,casei
2853  double precision, optional, intent(in) :: Ferr(0:rank)
2854  double complex, intent(in) :: MomVec(0:3,1:5), MomInv(15), masses2(0:5)
2855  double complex, intent(in) :: Ften(0:rank,0:rank,0:rank,0:rank)
2856  double complex, intent(in) :: Ftenuv(0:rank,0:rank,0:rank,0:rank)
2857  integer :: r,n0,n1,n2,n3,mode
2858  integer,parameter :: rankuv=8
2859  character(len=99) :: fname
2860  character(len=*),parameter :: fmt1 = "(A17,' = ',es23.16,' + i*',es23.16)"
2861  character(len=*),parameter :: fmt2 = "(A17,' = ',es23.16)"
2862  character(len=*),parameter :: fmt10 = "(' Ften(',i1,3(',',i1),') = ',es23.16,' + i*',es23.16)"
2863  character(len=*),parameter :: fmt11 = "(' Ferr(',i1,') = ',es23.16)"
2864  character(len=*),parameter :: fmt12 = "(' F0 = ',es23.16,' + i*',es23.16)"
2865 
2866 ! output file for result
2867  call getmode_cll(mode)
2868  select case (mode)
2869  case(1)
2870  fname = 'demo_6point_example00_coli.dat'
2871  case(2)
2872  fname = 'demo_6point_example00_dd.dat'
2873  case(3)
2874  fname = 'demo_6point_example00_comp.dat'
2875  end select
2876  write(fname(20:21),'(i2.2)') casei
2877 
2878 ! write(*,*) casei,fname
2879 
2880  open(unit=50,file=trim(fname),status='unknown')
2881 
2882  call getmuuv2_cll(muuv2)
2883  call getmuir2_cll(muir2)
2884  call getdeltauv_cll(deltauv)
2885  call getdeltair_cll(deltair1,deltair2)
2886 
2887  write (50,'(a37,i2,i3/)') ' Result for 6-point function, example',casen,casei
2888  write (50,'(a63,i2,i3,a)') ' The corresponding code can be found in demo.f90 under ''example',casen,casei,''''
2889  write (50,'((a))')
2890  write (50,'((a))') ' p42 '
2891  write (50,'((a))') ' ---------^--------- '
2892  write (50,'((a))') ' / \ '
2893  write (50,'((a))') ' m32,p3vec '
2894  write (50,'((a))') ' / p32 --------------------- p43 \ '
2895  write (50,'((a))') ' / / 3 \ \ '
2896  write (50,'((a))') ' p31 | m22,p2vec/2 4\ m42,p4vec | p53'
2897  write (50,'((a))') ' \ / \ / '
2898  write (50,'((a))') ' > p21 ----< >---- p54 < '
2899  write (50,'((a))') ' / \ / \ '
2900  write (50,'((a))') ' p20 | m12,p1vec\1 5/ m52,p5vec | p40'
2901  write (50,'((a))') ' \ \ 0 / / '
2902  write (50,'((a))') ' \ p10 --------------------- p50 / '
2903  write (50,'((a))') ' m02 '
2904  write (50,'((a))') ' \ / '
2905  write (50,'((a))') ' ---------v--------- '
2906  write (50,'((a))') ' p51 '
2907  write (50,'((a))') ''
2908  write (50,'((a))') ' Input:'
2909  write (50,fmt1) ' p1vec(0) ',momvec(0,1)
2910  write (50,fmt1) ' p1vec(1) ',momvec(1,1)
2911  write (50,fmt1) ' p1vec(2) ',momvec(2,1)
2912  write (50,fmt1) ' p1vec(3) ',momvec(3,1)
2913  write (50,fmt1) ' p2vec(0) ',momvec(0,2)
2914  write (50,fmt1) ' p2vec(1) ',momvec(1,2)
2915  write (50,fmt1) ' p2vec(2) ',momvec(2,2)
2916  write (50,fmt1) ' p2vec(3) ',momvec(3,2)
2917  write (50,fmt1) ' p3vec(0) ',momvec(0,3)
2918  write (50,fmt1) ' p3vec(1) ',momvec(1,3)
2919  write (50,fmt1) ' p3vec(2) ',momvec(2,3)
2920  write (50,fmt1) ' p3vec(3) ',momvec(3,3)
2921  write (50,fmt1) ' p4vec(0) ',momvec(0,4)
2922  write (50,fmt1) ' p4vec(1) ',momvec(1,4)
2923  write (50,fmt1) ' p4vec(2) ',momvec(2,4)
2924  write (50,fmt1) ' p4vec(3) ',momvec(3,4)
2925  write (50,fmt1) ' p5vec(0) ',momvec(0,5)
2926  write (50,fmt1) ' p5vec(1) ',momvec(1,5)
2927  write (50,fmt1) ' p5vec(2) ',momvec(2,5)
2928  write (50,fmt1) ' p5vec(3) ',momvec(3,5)
2929  write (50,fmt1) ' p10 ',mominv(1)
2930  write (50,fmt1) ' p21 ',mominv(2)
2931  write (50,fmt1) ' p32 ',mominv(3)
2932  write (50,fmt1) ' p43 ',mominv(4)
2933  write (50,fmt1) ' p54 ',mominv(5)
2934  write (50,fmt1) ' p50 ',mominv(6)
2935  write (50,fmt1) ' p20 ',mominv(7)
2936  write (50,fmt1) ' p31 ',mominv(8)
2937  write (50,fmt1) ' p42 ',mominv(9)
2938  write (50,fmt1) ' p53 ',mominv(10)
2939  write (50,fmt1) ' p40 ',mominv(11)
2940  write (50,fmt1) ' p51 ',mominv(12)
2941  write (50,fmt1) ' p30 ',mominv(13)
2942  write (50,fmt1) ' p41 ',mominv(14)
2943  write (50,fmt1) ' p52 ',mominv(15)
2944  write (50,fmt1) ' m02 ',masses2(0)
2945  write (50,fmt1) ' m12 ',masses2(1)
2946  write (50,fmt1) ' m22 ',masses2(2)
2947  write (50,fmt1) ' m32 ',masses2(3)
2948  write (50,fmt1) ' m42 ',masses2(4)
2949  write (50,fmt1) ' m52 ',masses2(5)
2950  if (rank.ge.rankuv) then
2951  write (50,fmt2) ' muUV2 ',muuv2
2952  end if
2953  write (50,fmt2) ' muIR2 ',muir2
2954  if (rank.ge.rankuv) then
2955  write (50,fmt2) ' deltaUV ',deltauv
2956  end if
2957  write (50,fmt2) ' DeltaIR1',deltair2
2958  write (50,fmt2) ' DeltaIR2',deltair1
2959  write (50,'((a))') ''
2960  write (50,'((a))') ' Conventions:'
2961  write (50,'((a))') ''
2962  if(rank.eq.0) then
2963  write (50,'((a))') ' (2*pi*mu)^(4-D) '
2964  write (50,'((a))') ' F0 = --------------- \int d^D q f(q,p_i)'
2965  write (50,'((a))') ' i*pi^2 '
2966 
2967  write (50,'((a))')
2968  write (50,'((a))') ' = F0_fin(muUV2,muIR2) + '// &
2969  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
2970  elseif(rank.ge.rankuv) then
2971  write (50,'((a))') ' (2*pi*mu)^(4-D) '
2972  write (50,'((a))') ' F = --------------- \int d^D q f(q,p_i)'
2973  write (50,'((a))') ' i*pi^2 '
2974 
2975  write (50,'((a))')
2976  write (50,'((a))') ' = F_fin(muUV2,muIR2) + a_UV*DeltaUV + '// &
2977  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
2978  else
2979  write (50,'((a))') ' (2*pi*mu)^(4-D) '
2980  write (50,'((a))') ' F = --------------- \int d^D q f(q,p_i)'
2981  write (50,'((a))') ' i*pi^2 '
2982 
2983  write (50,'((a))')
2984  write (50,'((a))') ' = F_fin(muUV2,muIR2) + '// &
2985  ' a_IR2*[DeltaIR2 + DeltaIR1*ln(muIR2)] + a_IR1*DeltaIR1'
2986  end if
2987  write (50,'((a))')
2988  write (50,'((a))') ' where'
2989  write (50,'((a))')
2990  if (rank.ge.rankuv) then
2991  write (50,'((a))') &
2992  ' c(epsUV) c(epsIR) c(epsIR)'
2993  write (50,'((a))') &
2994  ' DeltaUV = --------, DeltaIR1 = --------, DeltaIR2 = --------'
2995  write (50,'((a))') &
2996  ' epsUV epsIR epsIR^2'
2997  else
2998  write (50,'((a))') &
2999  ' c(epsIR) c(epsIR)'
3000  write (50,'((a))') &
3001  ' DeltaIR1 = --------, DeltaIR2 = --------'
3002  write (50,'((a))') &
3003  ' epsIR epsIR^2'
3004  end if
3005  write (50,'((a))')
3006  write (50,'((a))') ' c(eps) = (4*pi)^eps\Gamma(1+eps), D = 4 -2*eps '
3007  write (50,'((a))')
3008  write (50,'((a))') ' you can freely choose the regularization parameters'
3009  if (rank.ge.rankuv) then
3010  write (50,'((a))') ' of UV origin: muUV2 = mu^2, DeltaUV '
3011  end if
3012  write (50,'((a))') ' of IR origin: muIR2 = mu^2, DeltaIR1, DeltaIR2'
3013  write (50,'((a))')
3014  write (50,'((a))') ' note:'
3015  write (50,'((a))') ' - we effectively factor out a factor c(eps) '
3016  if (rank.ge.rankuv) then
3017  write (50,'((a))') ' - by default DeltaUV = DeltaIR1 = DeltaIR2 = 0 '
3018  else
3019  write (50,'((a))') ' - by default DeltaIR1 = DeltaIR2 = 0 '
3020  end if
3021  write (50,'((a))') ' - suitable DeltaIR2 can be used to adapt the effective normalization'
3022 
3023  write (50,'((a))') ''
3024  write (50,'((a))') ' Results:'
3025 
3026  if (rank.gt.0) then
3027  do r=0,rank
3028  do n0=0,r
3029  do n1=0,r-n0
3030  do n2=0,r-n0-n1
3031  n3=r-n0-n1-n2
3032  write (50,fmt10) n0,n1,n2,n3,ften(n0,n1,n2,n3)
3033  end do
3034  end do
3035  end do
3036  end do
3037 
3038  write (50,'(/(a))') ' Error estimates:'
3039  do r=0,rank
3040  write (50,fmt11) r,ferr(r)
3041  end do
3042  else
3043  write (50,fmt12) ften(0,0,0,0)
3044  end if
3045 
3046  write(*,'(/(a),(a)/)') ' The result has been written to the file ' &
3047  ,trim(fname)
3048 
3049  end subroutine writeresultften
3050 
3051 
3052 
3053 
3054 
3055 
3056  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3057  ! subroutine writeresultDB(caseN,casei,DBcoeff,DBcoeffuv,MomInv,masses2,rankm,rank,DBerr)
3058  !
3059  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3060 
3061  subroutine writeresultdb(caseN,casei,DBcoeff,DBcoeffuv,MomInv,masses2,rankm,rank,DBerr)
3062  implicit none
3063 
3064  integer, intent(in) :: rank,rankm,caseN,casei
3065  double complex, intent(in) :: MomInv(1), masses2(0:1)
3066  double complex, intent(in) :: DBcoeff(0:rank/2,0:rank)
3067  double complex, intent(in) :: DBcoeffuv(0:rank/2,0:rank)
3068  double precision, optional, intent(in) :: DBerr(0:rank)
3069  integer :: r,n0,n1,mode
3070  integer,parameter :: rankuv=0
3071  character(len=99) :: fname
3072  character(len=*),parameter :: fmt1 = "(A9,' = ',es23.16,' + i*',es23.16)"
3073  character(len=*),parameter :: fmt2 = "(A9,' = ',es23.16)"
3074  character(len=*),parameter :: fmt10 = "(' DBcoeff(',i1,1(',',i1),') = ',es23.16,' + i*',es23.16)"
3075  character(len=*),parameter :: fmt11 = "(' DBerr(',i1,') = ',es23.16)"
3076  character(len=*),parameter :: fmt12 = "(' DB0 = ',es23.16,' + i*',es23.16)"
3077  character(len=*),parameter :: fmt13 = "(' DB1 = ',es23.16,' + i*',es23.16)"
3078  character(len=*),parameter :: fmt14 = "(' DB00 = ',es23.16,' + i*',es23.16)"
3079  character(len=*),parameter :: fmt20 = "(' DBcoeffuv(',i1,1(',',i1),') = ',es23.16,' + i*',es23.16)"
3080 
3081 ! output file for result
3082  call getmode_cll(mode)
3083  select case (mode)
3084  case(1)
3085  fname = 'demo_2point_derivative_example00_coli.dat'
3086  case(2)
3087  fname = 'demo_2point_derivative_example00_dd.dat'
3088  case(3)
3089  fname = 'demo_2point_derivative_example00_comp.dat'
3090  end select
3091  write(fname(31:32),'(i2.2)') casei
3092 
3093  open(unit=50,file=trim(fname),status='unknown')
3094 
3095  call getmuuv2_cll(muuv2)
3096  call getmuir2_cll(muir2)
3097  call getdeltauv_cll(deltauv)
3098  call getdeltair_cll(deltair1,deltair2)
3099 
3100  write (50,'(a49,i2,i3/)') ' Result for 2-point function derivative, example ',casen,casei
3101  write (50,'((a))') ' m12 '
3102  write (50,'((a))') ' ------- '
3103  write (50,'((a))') ' / 1 \ '
3104  write (50,'((a))') ' / \ '
3105  write (50,'((a))') ' p10 ----- ----- p10 '
3106  write (50,'((a))') ' \ / '
3107  write (50,'((a))') ' \ 0 / '
3108  write (50,'((a))') ' ------- '
3109  write (50,'((a))') ' m02 '
3110  write (50,'((a))') ''
3111  write (50,'((a))') ' Input:'
3112  write (50,fmt1) ' p10 ',p10
3113  write (50,fmt1) ' m02 ',masses2(0)
3114  write (50,fmt1) ' m12 ',masses2(1)
3115  if (rank.ge.rankuv) then
3116  write (50,fmt2) ' muUV2 ',muuv2
3117  end if
3118 ! write (50,fmt2) ' muIR2 ',muIR2
3119  if (rank.ge.rankuv) then
3120  write (50,fmt2) ' deltaUV ',deltauv
3121  end if
3122 ! write (50,fmt2) ' deltaIR1',deltaIR2
3123 ! write (50,fmt2) ' deltaIR2',deltaIR1
3124  write (50,'((a))') ''
3125  write (50,'((a))') ' Conventions:'
3126  write (50,'((a))') ''
3127  if(rankm.eq.0.and.rank.gt.0) then
3128  write (50,'((a))') ' d (2*pi*mu)^(4-D) '
3129  write (50,'((a))') ' DB = ----- --------------- \int d^D q f(q,p_i)'
3130  write (50,'((a))') ' d p^2 i*pi^2 '
3131 
3132  write (50,'((a))')
3133  write (50,'((a))') ' = DB_fin(muUV2,muIR2) + a_UV*DeltaUV + a_IR1*DeltaIR1 '
3134  else if (rank.eq.0) then
3135  write (50,'((a))') ' d (2*pi*mu)^(4-D) '
3136  write (50,'((a))') ' DB0 = ----- --------------- \int d^D q f(q,p_i)'
3137  write (50,'((a))') ' d p^2 i*pi^2 '
3138 
3139  write (50,'((a))')
3140  write (50,'((a))') ' = DB0_fin(muIR2) + a_IR1*DeltaIR1 '
3141  else if (rank.eq.1.and.rankm.eq.1) then
3142  write (50,'((a))') ' d (2*pi*mu)^(4-D) '
3143  write (50,'((a))') ' DB1 = ----- --------------- \int d^D q f(q,p_i)'
3144  write (50,'((a))') ' d p^2 i*pi^2 '
3145 
3146  write (50,'((a))')
3147  write (50,'((a))') ' = DB1_fin(muIR2) + a_IR1*DeltaIR1 '
3148  else if (rank.eq.2.and.rankm.eq.2) then
3149  write (50,'((a))') ' d (2*pi*mu)^(4-D) '
3150  write (50,'((a))') ' DB00 = ----- --------------- \int d^D q f(q,p_i)'
3151  write (50,'((a))') ' d p^2 i*pi^2 '
3152 
3153  write (50,'((a))')
3154  write (50,'((a))') ' = DB0_fin(muUV2,muIR2)) + a_UV*DeltaUV + a_IR1*DeltaIR1 '
3155  end if
3156  write (50,'((a))')
3157  write (50,'((a))') ' where'
3158  write (50,'((a))')
3159  if (rank.ge.2) then
3160  write (50,'((a))') &
3161  ' c(epsUV) c(epsIR) '
3162  write (50,'((a))') &
3163  ' DeltaUV = --------, DeltaIR1 = -------- '
3164  write (50,'((a))') &
3165  ' epsUV epsIR '
3166  else
3167  write (50,'((a))') &
3168  ' c(epsIR) '
3169  write (50,'((a))') &
3170  ' DeltaIR1 = -------- '
3171  write (50,'((a))') &
3172  ' epsIR '
3173  end if
3174  write (50,'((a))')
3175  write (50,'((a))') ' c(eps) = (4*pi)^eps\Gamma(1+eps), D = 4 -2*eps '
3176  write (50,'((a))')
3177  write (50,'((a))') ' you can freely choose the regularization parameters'
3178  if (rank.ge.2) then
3179  write (50,'((a))') ' of UV origin: muUV2 = mu^2, DeltaUV '
3180  end if
3181  write (50,'((a))') ' of IR origin: muIR2 = mu^2, DeltaIR1'
3182  write (50,'((a))')
3183  write (50,'((a))') ' note:'
3184  write (50,'((a))') ' - we effectively factor out a factor c(eps) '
3185  write (50,'((a))') ' - by default DeltaUV = 0 '
3186 
3187  write (50,'((a))') ''
3188  write (50,'((a))') ' Results:'
3189 
3190  if (rank.gt.0) then
3191  if(rankm.eq.0) then
3192  do r=rankm,rank
3193  do n0=0,r/2
3194  n1=r-2*n0
3195  write (50,fmt10) n0,n1,dbcoeff(n0,n1)
3196  end do
3197  end do
3198  else if(rankm.eq.1.and.rank.eq.1) then
3199  write (50,fmt12) dbcoeff(0,1)
3200  else if(rankm.eq.2.and.rank.eq.2) then
3201  write (50,fmt13) dbcoeff(1,0)
3202  else
3203  write (*,*) 'writeresultDB: case not supported'
3204  end if
3205 
3206 
3207 ! do r=0,rank
3208 ! do n0=0,r/2
3209 ! n1=r-2*n0
3210 ! write (50,fmt13) n0,n1,DBcoeffuv(n0,n1)
3211 ! end do
3212 ! end do
3213 !
3214  if(present(dberr)) then
3215  write (50,'(/(a))') ' Error estimates:'
3216  do r=0,rank
3217  write (50,fmt11) r,dberr(r)
3218  end do
3219  end if
3220  else
3221  write (50,fmt12) dbcoeff(0,0)
3222  end if
3223 
3224  write(*,'(/(a),(a)/)') ' The result has been written to the file ' &
3225  ,trim(fname)
3226 
3227  end subroutine writeresultdb
3228 
3229 
3230  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3231  ! subroutine writeresultDB0(caseN,casei,DB0,MomInv,masses2)
3232  !
3233  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3234 
3235  subroutine writeresultdb0(caseN,casei,DB0,MomInv,masses2)
3236  implicit none
3237 
3238  integer, intent(in) :: caseN,casei
3239  double complex, intent(in) :: DB0
3240  double complex, intent(in) :: MomInv(1), masses2(0:1)
3241  double complex :: DBcoeff(0:0,0:0)
3242  double complex :: DBcoeffuv(0:0,0:0)
3243 
3244 
3245  dbcoeff(0:0,0:0) = db0
3246  dbcoeffuv(0:0,0:0) = 0d0
3247  call writeresultdb(casen,casei,dbcoeff,dbcoeffuv,mominv,masses2,0,0)
3248 
3249  end subroutine writeresultdb0
3250 
3251 
3252 
3253 
3254 
3255 
3256  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3257  ! subroutine getinvariants(N,MomInv,MomVec)
3258  !
3259  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3260 
3261  subroutine getinvariants(N,MomInv,MomVec)
3263  implicit none
3264 
3265  integer, intent(in) :: N
3266  double complex, intent(out) :: MomInv(N*(N-1)/2)
3267  double complex, intent(out), optional :: MomVec(0:3,N-1)
3268  double precision :: p(1:6,0:3)
3269 ! double complex :: MomInv6(1:15)
3270  double precision :: s(0:5,0:5)
3271  integer :: k,i
3272 
3273  if (n.gt.6.or.n.lt.3) then
3274  write(*,*) ' getinvariants: N>6 of N<3 not supported '
3275  else
3276 ! use momentum set from Racoon
3277  p(1,0) = -0.9500000000000000d+02
3278  p(1,1) = 0.0000000000000000d+00
3279  p(1,2) = 0.0000000000000000d+00
3280  p(1,3) = 0.9500000000000000d+02
3281 
3282  p(2,0) = -0.9500000000000000d+02
3283  p(2,1) = 0.0000000000000000d+00
3284  p(2,2) = 0.0000000000000000d+00
3285  p(2,3) = -0.9500000000000000d+02
3286 
3287  p(3,0) = 0.2046111001757171d+02
3288  p(3,1) = 0.1057734233089455d+02
3289  p(3,2) = -0.2324961261504543d+01
3290  p(3,3) = 0.1736005205921753d+02
3291 
3292  p(4,0) = 0.3558305163378869d+01
3293  p(4,1) = 0.1436222934374051d+01
3294  p(4,2) = -0.2174258125294355d+01
3295  p(4,3) = -0.2423097382091398d+01
3296 
3297  p(5,0) = 0.8154540918019539d+02
3298  p(5,1) = -0.5230395944682889d+02
3299  p(5,2) = 0.3083642435466509d+02
3300  p(5,3) = 0.5443403822581044d+02
3301 
3302  p(6,0) = 0.8443517563885433d+02
3303  p(6,1) = 0.4029039418156027d+02
3304  p(6,2) = -0.2633720496786619d+02
3305  p(6,3) = -0.6937099290293661d+02
3306 
3307 
3308  do i=n+1,6
3309  p(n,0) = p(n,0) + p(i,0)
3310  p(n,1) = p(n,1) + p(i,1)
3311  p(n,2) = p(n,2) + p(i,2)
3312  p(n,3) = p(n,3) + p(i,3)
3313  end do
3314 
3315  if (present(momvec)) then
3316  momvec(0:3,1) = p(1,0:3)
3317  do i=2,n-1
3318  momvec(0:3,i) = momvec(0:3,i-1)+p(i,0:3)
3319  end do
3320  end if
3321 
3322  ! define invariants
3323  do k=0,n-1
3324  s(modulo(k+1,n),k) = p(k+1,0)**2
3325  do i=1,3
3326  s(modulo(k+1,n),k) = s(modulo(k+1,n),k) - p(k+1,i)**2
3327  end do
3328  if(abs(s(modulo(k+1,n),k)).lt.1d-14* abs(p(k+1,0))**2) s(modulo(k+1,n),k) = 0d0
3329 
3330  s(k,modulo(k+1,n)) = s(modulo(k+1,n),k)
3331  end do
3332  do k=0,n-1
3333  s(modulo(k+2,n),k) = (p(modulo(k+1,n)+1,0) + p(k+1,0))**2
3334  do i=1,3
3335  s(modulo(k+2,n),k) = s(modulo(k+2,n),k) &
3336  - (p(modulo(k+1,n)+1,i) + p(k+1,i))**2
3337  end do
3338  s(k,modulo(k+2,n)) = s(modulo(k+2,n),k)
3339  end do
3340  do k=0,(n-1)/2
3341  s(modulo(k+3,n),k) = (p(modulo(k+2,n)+1,0) + p(modulo(k+1,n)+1,0) + p(k+1,0))**2
3342  do i=1,3
3343  s(modulo(k+3,n),k) = s(modulo(k+3,n),k) &
3344  - (p(modulo(k+2,n)+1,i) + p(modulo(k+1,n)+1,i) + p(k+1,i))**2
3345  end do
3346  s(k,modulo(k+3,n)) = s(modulo(k+3,n),k)
3347  end do
3348 
3349  do i=0,(n+1)/2-2
3350  do k=1,n
3351  mominv(k+n*i) = s(modulo(k+i,n),k-1)
3352 
3353 ! write(*,*) 'mominv',k+N*i,s(modulo(k+i,N),k-1)
3354  end do
3355  end do
3356  if(modulo(n,2).eq.0) then
3357  do k=1,n/2
3358  mominv(k+n*(n/2-1)) = s(modulo(k+n/2-1,n),k-1)
3359 ! write(*,*) 'mominv',k+N*i,s(modulo(k+N/2-1,N),k-1)
3360  end do
3361  end if
3362  end if
3363 
3364  end subroutine getinvariants
3365 
3366 
3367 end program demo
writeresulteten
subroutine writeresulteten(caseN, casei, Eten, Etenuv, MomVec, MomInv, masses2, rank, Eerr)
Definition: demo.f90:2430
writeresultb
subroutine writeresultb(caseN, casei, Bcoeff, Bcoeffuv, MomInv, masses2, rank, Berr)
Definition: demo.f90:1319
writeresultd0
subroutine writeresultd0(caseN, casei, D0, MomInv, masses2)
Definition: demo.f90:2012
writeresultc
subroutine writeresultc(caseN, casei, Ccoeff, Ccoeffuv, MomInv, masses2, rank, Cerr)
Definition: demo.f90:1477
writeresulte0
subroutine writeresulte0(caseN, casei, E0, MomInv, masses2)
Definition: demo.f90:2406
writeresultf0
subroutine writeresultf0(caseN, casei, F0, MomInv, masses2)
Definition: demo.f90:2824
collier
Definition: COLLIER.F90:30
demo
program demo
Definition: demo.f90:10
writeresultd
subroutine writeresultd(caseN, casei, Dcoeff, Dcoeffuv, MomInv, masses2, rank, Derr)
Definition: demo.f90:1840
writeresultften
subroutine writeresultften(caseN, casei, Ften, Ftenuv, MomVec, MomInv, masses2, rank, Ferr)
Definition: demo.f90:2850
writeresultdten
subroutine writeresultdten(caseN, casei, Dten, Dtenuv, MomVec, MomInv, masses2, rank, Derr)
Definition: demo.f90:2036
writeresultf
subroutine writeresultf(caseN, casei, Fcoeff, Fcoeffuv, MomInv, masses2, rank, Ferr)
Definition: demo.f90:2632
writeresulta
subroutine writeresulta(caseN, casei, Acoeff, Acoeffuv, masses2, rank, Aerr)
Definition: demo.f90:1179
writeresulte
subroutine writeresulte(caseN, casei, Ecoeff, Ecoeffuv, MomInv, masses2, rank, Eerr)
Definition: demo.f90:2222
writeresultc0
subroutine writeresultc0(caseN, casei, C0, MomInv, masses2)
Definition: demo.f90:1641
writeresultcten
subroutine writeresultcten(caseN, casei, Cten, Ctenuv, MomVec, MomInv, masses2, rank, Cerr)
Definition: demo.f90:1664
writeresultb0
subroutine writeresultb0(caseN, casei, B0, MomInv, masses2)
Definition: demo.f90:1454
writeresultdb
subroutine writeresultdb(caseN, casei, DBcoeff, DBcoeffuv, MomInv, masses2, rankm, rank, DBerr)
Definition: demo.f90:3062
writeresulta0
subroutine writeresulta0(caseN, casei, A0, masses2)
Definition: demo.f90:1297
writeresultdb0
subroutine writeresultdb0(caseN, casei, DB0, MomInv, masses2)
Definition: demo.f90:3236
getinvariants
subroutine getinvariants(N, MomInv, MomVec)
Definition: demo.f90:3262