-
Notifications
You must be signed in to change notification settings - Fork 9
/
MSbarmass.f
504 lines (502 loc) · 13.4 KB
/
MSbarmass.f
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
************************************************************************
*
* MSbarmass.f:
*
* Function that provides the running of the heavy quark masses in
* the MSbar scheme.
*
************************************************************************
FUNCTION MSBARMASS(IM,Q2)
*
IMPLICIT NONE
*
include "../commons/m2th.h"
include "../commons/kren.h"
include "../commons/Evs.h"
include "../commons/MaxFlavourAlpha.h"
include "../commons/Nf_FF.h"
include "../commons/ThresholdAlphaQCD.h"
**
* Input Variables
*
INTEGER IM
DOUBLE PRECISION Q2
**
* Internal Variables
*
DOUBLE PRECISION A_QCD,ASQ,ASI,ASIM,ASTH(4:6),ASTHM(4:6)
DOUBLE PRECISION INMASS,EVF
DOUBLE PRECISION EVMASS,DECOUP
DOUBLE PRECISION EPS
PARAMETER(EPS=1D-10)
**
* Output Variables
*
DOUBLE PRECISION MSBARMASS
*
IF(IM.LT.4.OR.IM.GT.6)THEN
WRITE(6,*) "In src/Evolution/MSbarmass.f:"
WRITE(6,*) "IM out of range, IM =",IM
CALL EXIT(-10)
ENDIF
*
ASQ = A_QCD(Q2)
*
ASTH(4) = asthUp(4)
ASTH(5) = asthUp(5)
ASTH(6) = asthUp(6)
*
ASTHM(4) = asthDown(4)
ASTHM(5) = asthDown(5)
ASTHM(6) = asthDown(6)
*
ASI = asthUp(IM)
ASIM = asthDown(IM)
*
INMASS = DSQRT(M2PH(IM))
*
IF(EVS.EQ."FF")THEN
EVF = EVMASS(NF_FF,ASI,ASQ)
ELSEIF(EVS.EQ."VF")THEN
* Charm
IF(IM.EQ.4)THEN
IF(Q2.GE.M2TH(6).AND.NFMAXALPHA.GE.6)THEN
EVF = EVMASS(4,ASI,ASTHM(5))
1 * DECOUP("UP",5,DLOG(K2TH(5)))
2 * EVMASS(5,ASTH(5),ASTHM(6))
3 * DECOUP("UP",6,DLOG(K2TH(6)))
4 * EVMASS(6,ASTH(6),ASQ)
ELSEIF(Q2.GE.M2TH(5).AND.NFMAXALPHA.GE.5)THEN
EVF = EVMASS(4,ASI,ASTHM(5))
1 * DECOUP("UP",5,DLOG(K2TH(5)))
2 * EVMASS(5,ASTH(5),ASQ)
ELSEIF(Q2.GE.M2TH(4).AND.NFMAXALPHA.GE.4)THEN
EVF = EVMASS(4,ASI,ASQ)
ELSE
EVF = DECOUP("DW",4,DLOG(K2TH(4))) / EVMASS(3,ASQ,ASIM)
ENDIF
* Bottom
ELSEIF(IM.EQ.5)THEN
IF(Q2.GE.M2TH(6).AND.NFMAXALPHA.GE.6)THEN
EVF = EVMASS(5,ASI,ASTHM(6))
1 * DECOUP("UP",6,DLOG(K2TH(6)))
2 * EVMASS(6,ASTH(6),ASQ)
ELSEIF(Q2.GE.M2TH(5).AND.NFMAXALPHA.GE.5)THEN
EVF = EVMASS(5,ASI,ASQ)
ELSEIF(Q2.GE.M2TH(4).AND.NFMAXALPHA.GE.4)THEN
EVF = DECOUP("DW",5,DLOG(K2TH(5))) / EVMASS(4,ASQ,ASIM)
ELSE
EVF = DECOUP("DW",4,DLOG(K2TH(4)))
1 / EVMASS(3,ASQ,ASTHM(4))
2 * DECOUP("DW",5,DLOG(K2TH(5)))
3 / EVMASS(4,ASTH(4),ASIM)
ENDIF
* Top
ELSEIF(IM.EQ.6)THEN
IF(Q2.GE.M2TH(6).AND.NFMAXALPHA.GE.6)THEN
EVF = EVMASS(6,ASI,ASQ)
ELSEIF(Q2.GE.M2TH(5).AND.NFMAXALPHA.GE.5)THEN
EVF = DECOUP("DW",6,DLOG(K2TH(6))) / EVMASS(5,ASQ,ASIM)
ELSEIF(Q2.GE.M2TH(4).AND.NFMAXALPHA.GE.4)THEN
EVF = DECOUP("DW",5,DLOG(K2TH(5)))
1 / EVMASS(4,ASQ,ASTHM(5))
2 * DECOUP("DW",6,DLOG(K2TH(6)))
3 / EVMASS(5,ASTH(5),ASIM)
ELSE
EVF = DECOUP("DW",4,DLOG(K2TH(4)))
1 / EVMASS(3,ASQ,ASTHM(4))
2 * DECOUP("DW",5,DLOG(K2TH(5)))
3 / EVMASS(4,ASTH(4),ASTHM(5))
4 * DECOUP("DW",6,DLOG(K2TH(6)))
5 / EVMASS(5,ASTH(5),ASIM)
ENDIF
ENDIF
ENDIF
*
MSBARMASS = INMASS * EVF
*
RETURN
END
*
****************************************************************
FUNCTION EVMASS(NF,AS0,AS)
*
IMPLICIT NONE
*
include "../commons/ipt.h"
include "../commons/AlphaEvolution.h"
**
* Input Variables
*
INTEGER NF
DOUBLE PRECISION AS,AS0
**
* Internal Variables
*
DOUBLE PRECISION FACT,BT0
DOUBLE PRECISION BETA0APF,BETA1APF,BETA2APF,B1,B2
DOUBLE PRECISION GAMMA0APF,GAMMA1APF,GAMMA2APF,C0,C1,C2
DOUBLE PRECISION INTEGRANDMASSRUNNING,DGAUSS,EPS
PARAMETER(EPS=1D-7)
EXTERNAL INTEGRANDMASSRUNNING
INTEGER CMNF,CMIPT
COMMON / CMASSRUNNING / CMNF,CMIPT
**
* Output Variables
*
DOUBLE PRECISION EVMASS
*
IF(AlphaEvol(1:8).EQ."expanded")THEN
BT0 = BETA0APF(NF)
*
C0 = GAMMA0APF() / BT0
FACT = DEXP( C0 * DLOG( AS / AS0 ) )
IF(IPT.EQ.0)THEN
EVMASS = FACT
ELSEIF(IPT.EQ.1)THEN
B1 = BETA1APF(NF) / BT0
C1 = GAMMA1APF(NF) / BT0
EVMASS = FACT * ( 1D0 + ( C1 - B1 * C0 ) * AS )
1 / ( 1D0 + ( C1 - B1 * C0 ) * AS0 )
ELSEIF(IPT.EQ.2)THEN
B1 = BETA1APF(NF) / BT0
B2 = BETA2APF(NF) / BT0
C1 = GAMMA1APF(NF) / BT0
C2 = GAMMA2APF(NF) / BT0
EVMASS = FACT * ( 1D0 + ( C1 - B1 * C0 ) * AS
1 + ( C2 - C1 * B1 - B2 * C0
2 + B1**2 * C0 + ( C1
3 - B1 * C0 )**2 ) * AS**2 / 2D0 )
4 / ( 1D0 + ( C1 - B1 * C0 ) * AS0
5 + ( C2 - C1 * B1 - B2 * C0
6 + B1**2 * C0 + ( C1
7 - B1 * C0 )**2 ) * AS0**2 / 2D0 )
ENDIF
ELSE
CMNF = NF
CMIPT = IPT
EVMASS = DEXP(DGAUSS(INTEGRANDMASSRUNNING,AS0,AS,EPS))
ENDIF
*
RETURN
END
*
****************************************************************
FUNCTION INTEGRANDMASSRUNNING(AS)
**
* Input Variables
*
DOUBLE PRECISION AS
**
* Internal Variables
*
DOUBLE PRECISION FGAMMA,FBETA
INTEGER CMNF,CMIPT
COMMON / CMASSRUNNING / CMNF,CMIPT
**
* Output Variables
*
DOUBLE PRECISION INTEGRANDMASSRUNNING
*
INTEGRANDMASSRUNNING = FGAMMA(AS,CMNF,CMIPT)
1 / FBETA(AS,CMNF,CMIPT)
*
RETURN
END
*
****************************************************************
FUNCTION DECOUP(DIR,IM,LN)
*
IMPLICIT NONE
*
include "../commons/ipt.h"
include "../commons/m2th.h"
include "../commons/kren.h"
include "../commons/ThresholdAlphaQCD.h"
**
* Input Variables
*
INTEGER IM
DOUBLE PRECISION LN
CHARACTER*2 DIR
**
* Internal Variables
*
DOUBLE PRECISION ASTH2
DOUBLE PRECISION EPS
PARAMETER(EPS=1D-7)
**
* Output Variables
*
DOUBLE PRECISION DECOUP
*
IF(IPT.LE.1)THEN
DECOUP = 1D0
RETURN
ELSE
IF(DIR.EQ."DW")THEN
ASTH2 = asthUp(IM)**2
DECOUP = 1D0 + ASTH2 * ( 89D0 / 27D0 - 20D0 / 9D0 * LN
1 + 4D0 / 3D0 * LN**2 )
ELSEIF(DIR.EQ."UP")THEN
ASTH2 = asthDown(IM)**2
DECOUP = 1D0 - ASTH2 * ( 89D0 / 27D0 - 20D0 / 9D0 * LN
1 + 4D0 / 3D0 * LN**2 )
ELSE
WRITE(6,*) "In src/Evolution/MSbarmass.f:"
WRITE(6,*) "Unknown direction, DIR =",DIR
CALL EXIT(-10)
ENDIF
ENDIF
*
RETURN
END
*
****************************************************************************
*
* QCD gamma function.
*
****************************************************************************
function fgamma(a,nf,ipt)
*
implicit none
**
* Input Variables
*
double precision a
integer nf,ipt
**
* Internal Variables
*
double precision gamma0apf,gamma1apf,gamma2apf,gamma3apf
**
* Output Variables
*
double precision fgamma
*
if(ipt.eq.0)then
fgamma = - a * gamma0apf()
elseif(ipt.eq.1)then
fgamma = - a * ( gamma0apf() + a * gamma1apf(nf) )
elseif(ipt.eq.2)then
fgamma = - a * ( gamma0apf()
1 + a * ( gamma1apf(nf) + a * gamma2apf(nf) ) )
elseif(ipt.eq.3)then
fgamma = - a * ( gamma0apf()
1 + a * ( gamma1apf(nf)
2 + a * ( gamma2apf(nf) + a * gamma3apf(nf) ) ) )
endif
*
return
end
*
****************************************************************************
function gamma0apf()
*
implicit none
**
* Output Variables
*
double precision gamma0apf
*
gamma0apf = 4d0
*
return
end
*
****************************************************************************
function gamma1apf(nf)
*
implicit none
**
* Input Variables
*
integer nf
**
* Output Variables
*
double precision gamma1apf
*
gamma1apf = 202d0 / 3d0 - 20d0 / 9d0 * nf
*
return
end
*
****************************************************************************
function gamma2apf(nf)
*
implicit none
*
include "../commons/consts.h"
**
* Input Variables
*
integer nf
**
* Output Variables
*
double precision gamma2apf
*
gamma2apf = 1249d0 - ( 2216d0 / 27d0 + 160d0 / 3d0 * zeta3 ) * nf
1 - 140d0 / 81d0 * nf**2
*
return
end
*
****************************************************************************
function gamma3apf(nf)
*
implicit none
*
include "../commons/consts.h"
**
* Input Variables
*
integer nf
**
* Output Variables
*
double precision gamma3apf
*
gamma3apf = 4603055d0 / 162d0 + 135680d0 * zeta3 / 28d0
1 - 8800d0 * zeta5
2 + ( - 91723d0 / 27d0 - 34192d0 * zeta3 / 9d0 + 880d0 * zeta4
3 + 18400d0 * zeta5 / 9d0 ) * nf
5 + ( 5242d0 / 243d0 + 800d0 * zeta3 / 9d0
6 - 160d0 * zeta4 / 3d0 ) * nf**2
7 + ( 332d0 / 243d0 + 64d0 * zeta3 / 27d0 ) * nf**3
*
return
end
*
****************************************************************************
*
* Function whose zero is the so-called RG invariant mass.
* For the MSbar mass, it finds the scale m such that m(m) = m.
*
****************************************************************************
function MassQSplit(i,Q)
*
implicit none
*
include "../commons/m2th.h"
include "../commons/kren.h"
include "../commons/Evs.h"
include "../commons/MaxFlavourAlpha.h"
include "../commons/Nf_FF.h"
include "../commons/ThresholdAlphaQCD.h"
**
* Input Variables
*
integer i
double precision Q
**
* Internal Variables
*
DOUBLE PRECISION Q2
DOUBLE PRECISION A_QCD,ASQ,ASI,ASTH(4:6),ASTHM(4:6)
DOUBLE PRECISION LN
DOUBLE PRECISION INMASS,EVF
DOUBLE PRECISION EVMASS,DECOUP
*
* Output Variables
*
double precision MassQSplit
*
Q2 = Q * Q
*
ASI = A_QCD(Q2)
ASQ = A_QCD(Q2TH(I))
*
ASTH(4) = asthUp(4)
ASTH(5) = asthUp(5)
ASTH(6) = asthUp(6)
*
ASTHM(4) = asthDown(4)
ASTHM(5) = asthDown(5)
ASTHM(6) = asthDown(6)
*
INMASS = DSQRT(M2Q(I))
*
IF(EVS.EQ."FF")THEN
EVF = EVMASS(NF_FF,ASI,ASQ)
ELSEIF(EVS.EQ."VF")THEN
LN = DLOG(KREN)
* Charm
IF(I.EQ.4)THEN
IF(Q2.GE.M2TH(6).AND.NFMAXALPHA.GE.6)THEN
EVF = EVMASS(6,ASI,ASQ)
ELSEIF(Q2.GE.M2TH(5).AND.NFMAXALPHA.GE.5)THEN
IF(Q2TH(4).GE.M2TH(6).AND.NFMAXALPHA.GE.6)THEN
EVF = EVMASS(5,ASI, ASTHM(6)) * DECOUP("UP",6,LN)
1 * EVMASS(6,ASTH(6),ASQ)
ELSE
EVF = EVMASS(5,ASI,ASQ)
ENDIF
ELSEIF(Q2.GE.M2TH(4).AND.NFMAXALPHA.GE.4)THEN
IF(Q2TH(4).GE.M2TH(6).AND.NFMAXALPHA.GE.6)THEN
EVF = EVMASS(4,ASI, ASTHM(5)) * DECOUP("UP",5,LN)
1 * EVMASS(5,ASTH(5), ASTHM(6)) * DECOUP("UP",6,LN)
2 * EVMASS(6,ASTH(6),ASQ)
ELSEIF(Q2TH(4).GE.M2TH(5).AND.NFMAXALPHA.GE.5)THEN
EVF = EVMASS(4,ASI, ASTHM(5)) * DECOUP("UP",5,LN)
2 * EVMASS(5,ASTH(5),ASQ)
ELSE
EVF = EVMASS(4,ASI,ASQ)
ENDIF
ENDIF
* Bottom
ELSEIF(I.EQ.5)THEN
IF(Q2.GE.M2TH(6).AND.NFMAXALPHA.GE.6)THEN
EVF = EVMASS(6,ASI,ASQ)
ELSEIF(Q2.GE.M2TH(5).AND.NFMAXALPHA.GE.5)THEN
IF(Q2TH(5).GE.M2TH(6).AND.NFMAXALPHA.GE.6)THEN
EVF = EVMASS(5,ASI, ASTHM(6)) * DECOUP("UP",6,LN)
1 * EVMASS(6,ASTH(6),ASQ)
ELSE
EVF = EVMASS(5,ASI,ASQ)
ENDIF
ENDIF
* Top
ELSEIF(I.EQ.6)THEN
EVF = EVMASS(6,ASI,ASQ)
ENDIF
ENDIF
*
MassQSplit = INMASS / EVF - Q
*
return
end
*
****************************************************************************
*
* Subroutine that computes the RG invariant masses
*
****************************************************************************
subroutine ComputeRGInvariantMasses
*
implicit none
*
include "../commons/m2th.h"
**
* Internal Variables
*
integer i
double precision MassQSplit,zriddr
double precision x1,x2
double precision acc
parameter(acc=1d-10)
external MassQSplit
*
do i=6,4,-1
call ThresholdAlphaQCD
if(q2th(i).ne.m2q(i))then
x1 = dsqrt(m2q(i))
x2 = dsqrt(q2th(i))
m2ph(i) = zriddr(MassQSplit,i,x1,x2,acc)**2
endif
call ComputeHeavyQuarkThresholds
enddo
*
return
end