Source code for ekore.anomalous_dimensions.unpolarized.space_like.as4.gqg

# pylint: skip-file
# fmt: off
r"""The unpolarized, space-like anomalous dimension
:math:`\gamma_{qg}^{(3)}`."""
import numba as nb
import numpy as np

from .....harmonics import cache as c
from .....harmonics.log_functions import (
    lm11,
    lm11m1,
    lm12,
    lm12m1,
    lm13,
    lm13m1,
    lm14,
    lm14m1,
    lm15,
    lm15m1,
)


[docs] @nb.njit(cache=True) def gamma_qg_nf3(n, cache): r"""Return the part proportional to :math:`nf^3` of :math:`\gamma_{qg}^{(3)}`. The expression is copied exact from :eqref:`3.12` of :cite:`Davies:2016jie`. Parameters ---------- n : complex Mellin moment cache: numpy.ndarray Harmonic sum cache Returns ------- complex |N3LO| non-singlet anomalous dimension :math:`\gamma_{qg}^{(3)}|_{nf^3}` """ S1 = c.get(c.S1, cache, n) S2 = c.get(c.S2, cache, n) Sm2 = c.get(c.Sm2, cache, n, is_singlet=True) S3 = c.get(c.S3, cache, n) S21 = c.get(c.S21, cache, n) Sm3 = c.get(c.Sm3, cache, n, is_singlet=True) S4 = c.get(c.S4, cache, n) Sm4 = c.get(c.Sm4, cache, n, is_singlet=True) S31 = c.get(c.S31, cache, n) S211 = c.get(c.S211, cache, n) return 1.3333333333333333 * ( 44.56685134331718 / (-1.0 + n) - 82.37037037037037 / np.power(n, 5) + 95.30864197530865 / np.power(n, 4) - 298.6951686088834 / np.power(n, 3) + 334.4519003852186 / np.power(n, 2) - 576.1641081960868 / n - 156.44444444444446 / np.power(1.0 + n, 6) + 271.4074074074074 / np.power(1.0 + n, 5) - 142.6172839506173 / np.power(1.0 + n, 4) + 49.20926725891911 / np.power(1.0 + n, 3) + 242.3162373306182 / np.power(1.0 + n, 2) + 383.95514040473176 / (1.0 + n) + 75.85185185185185 / np.power(2.0 + n, 5) - 18.962962962962962 / np.power(2.0 + n, 4) - 28.97119341563786 / np.power(2.0 + n, 3) + 57.904384241653375 / np.power(2.0 + n, 2) + 153.56036440484917 / (2.0 + n) - (7.725651577503429 * S1) / (-1.0 + n) + (35.55555555555556 * S1) / np.power(n, 5) - (53.333333333333336 * S1) / np.power(n, 4) + (149.7283950617284 * S1) / np.power(n, 3) - (189.49794238683128 * S1) / np.power(n, 2) + (219.77429764880566 * S1) / n - (71.11111111111111 * S1) / np.power(1.0 + n, 5) + (75.85185185185185 * S1) / np.power(1.0 + n, 4) + (45.4320987654321 * S1) / np.power(1.0 + n, 3) - (24.691358024691358 * S1) / np.power(1.0 + n, 2) - (242.01773110008048 * S1) / (1.0 + n) + (37.925925925925924 * S1) / np.power(2.0 + n, 4) + (53.72839506172839 * S1) / np.power(2.0 + n, 3) + (39.76954732510288 * S1) / np.power(2.0 + n, 2) + (33.84214810968268 * S1) / (2.0 + n) - (8.954732510288066 * S2) / (-1.0 + n) + (28.444444444444443 * S2) / np.power(n, 4) - (47.407407407407405 * S2) / np.power(n, 3) + (117.33333333333333 * S2) / np.power(n, 2) - (128.52674897119343 * S2) / n - (71.11111111111111 * S2) / np.power(1.0 + n, 4) + (9.481481481481481 * S2) / np.power(1.0 + n, 3) + (60.44444444444444 * S2) / np.power(1.0 + n, 2) - (11.588477366255145 * S2) / (1.0 + n) + (56.888888888888886 * S2) / np.power(2.0 + n, 3) + (112.19753086419753 * S2) / np.power(2.0 + n, 2) + (144.98765432098764 * S2) / (2.0 + n) - (2.3703703703703702 * (np.power(S1, 2) + S2)) / (-1.0 + n) + (2.3703703703703702 * (np.power(S1, 2) + S2)) / np.power(n, 3) + (3.3580246913580245 * (np.power(S1, 2) + S2)) / np.power(n, 2) - (13.695473251028806 * (np.power(S1, 2) + S2)) / n + (7.111111111111111 * (np.power(S1, 2) + S2)) / np.power(1.0 + n, 4) + (9.481481481481481 * (np.power(S1, 2) + S2)) / np.power(1.0 + n, 3) - (2.962962962962963 * (np.power(S1, 2) + S2)) / np.power(1.0 + n, 2) + (53.76131687242798 * (np.power(S1, 2) + S2)) / (1.0 + n) - (9.481481481481481 * (np.power(S1, 2) + S2)) / np.power(2.0 + n, 3) - (21.333333333333332 * (np.power(S1, 2) + S2)) / np.power(2.0 + n, 2) - (38.650205761316876 * (np.power(S1, 2) + S2)) / (2.0 + n) - (3.1604938271604937 * S21) / (-1.0 + n) + (7.111111111111111 * S21) / np.power(n, 3) - (10.666666666666666 * S21) / np.power(n, 2) + (32.0 * S21) / n - (14.222222222222221 * S21) / np.power(1.0 + n, 3) - (53.333333333333336 * S21) / (1.0 + n) + (9.481481481481481 * S21) / np.power(2.0 + n, 2) + (24.493827160493826 * S21) / (2.0 + n) - (3.1604938271604937 * S3) / (-1.0 + n) + (7.111111111111111 * S3) / np.power(n, 3) - (10.666666666666666 * S3) / np.power(n, 2) + (28.049382716049383 * S3) / n - (14.222222222222221 * S3) / np.power(1.0 + n, 3) - (43.06172839506173 * S3) / (1.0 + n) + (9.481481481481481 * S3) / np.power(2.0 + n, 2) + (14.222222222222221 * S3) / (2.0 + n) - (3.1604938271604937 * (S1 * S2 - 1.0 * S21 + S3)) / (-1.0 + n) + (7.111111111111111 * (S1 * S2 - 1.0 * S21 + S3)) / np.power(n, 3) - (10.666666666666666 * (S1 * S2 - 1.0 * S21 + S3)) / np.power(n, 2) + (32.0 * (S1 * S2 - 1.0 * S21 + S3)) / n - (14.222222222222221 * (S1 * S2 - 1.0 * S21 + S3)) / np.power(1.0 + n, 3) - (53.333333333333336 * (S1 * S2 - 1.0 * S21 + S3)) / (1.0 + n) + (9.481481481481481 * (S1 * S2 - 1.0 * S21 + S3)) / np.power(2.0 + n, 2) + (24.493827160493826 * (S1 * S2 - 1.0 * S21 + S3)) / (2.0 + n) + (0.5267489711934157 * (np.power(S1, 3) + 3.0 * S1 * S2 + 2.0 * S3)) / (-1.0 + n) - (1.1851851851851851 * (np.power(S1, 3) + 3.0 * S1 * S2 + 2.0 * S3)) / np.power(n, 3) + (1.9753086419753085 * (np.power(S1, 3) + 3.0 * S1 * S2 + 2.0 * S3)) / np.power(n, 2) - (4.674897119341564 * (np.power(S1, 3) + 3.0 * S1 * S2 + 2.0 * S3)) / n + (2.3703703703703702 * (np.power(S1, 3) + 3.0 * S1 * S2 + 2.0 * S3)) / np.power(1.0 + n, 3) + (7.57201646090535 * (np.power(S1, 3) + 3.0 * S1 * S2 + 2.0 * S3)) / (1.0 + n) - (1.5802469135802468 * (np.power(S1, 3) + 3.0 * S1 * S2 + 2.0 * S3)) / np.power(2.0 + n, 2) - (2.765432098765432 * (np.power(S1, 3) + 3.0 * S1 * S2 + 2.0 * S3)) / (2.0 + n) - ( 1.1851851851851851 * ( 0.041666666666666664 * np.power(S1, 4) + 0.25 * np.power(S1, 2) * S2 + 0.125 * np.power(S2, 2) + 0.3333333333333333 * S1 * S3 + 0.25 * S4 ) ) / n + ( 2.3703703703703702 * ( 0.041666666666666664 * np.power(S1, 4) + 0.25 * np.power(S1, 2) * S2 + 0.125 * np.power(S2, 2) + 0.3333333333333333 * S1 * S3 + 0.25 * S4 ) ) / (1.0 + n) - ( 2.3703703703703702 * ( 0.041666666666666664 * np.power(S1, 4) + 0.25 * np.power(S1, 2) * S2 + 0.125 * np.power(S2, 2) + 0.3333333333333333 * S1 * S3 + 0.25 * S4 ) ) / (2.0 + n) + (3.5555555555555554 * S4) / n - (7.111111111111111 * S4) / (1.0 + n) + (7.111111111111111 * S4) / (2.0 + n) ) + 3.0 * ( 2.5381463063368415 / (-1.0 + n) + 3.5555555555555554 / np.power(n, 5) + 5.728395061728395 / np.power(n, 4) + 4.559670781893004 / np.power(n, 3) + 2.036401939671693 / np.power(n, 2) - 1.1169664019346115 / n + 26.074074074074073 / np.power(1.0 + n, 5) - 56.098765432098766 / np.power(1.0 + n, 4) + 70.0246913580247 / np.power(1.0 + n, 3) - 56.58563233464838 / np.power(1.0 + n, 2) + 23.268759906664844 / (1.0 + n) + 18.962962962962962 / np.power(2.0 + n, 5) - 37.925925925925924 / np.power(2.0 + n, 4) + 22.386831275720166 / np.power(2.0 + n, 3) + 25.899964373170548 / np.power(2.0 + n, 2) - 33.71051594275432 / (2.0 + n) + (0.3511659807956104 * S1) / (-1.0 + n) - (9.481481481481481 * S1) / np.power(n, 4) + (8.938271604938272 * S1) / np.power(n, 3) - (17.432098765432098 * S1) / np.power(n, 2) + (15.015825807984449 * S1) / n + (14.222222222222221 * S1) / np.power(1.0 + n, 4) - (1.7777777777777777 * S1) / np.power(1.0 + n, 3) - (14.534979423868313 * S1) / np.power(1.0 + n, 2) + (48.28933603835209 * S1) / (1.0 + n) + (9.481481481481481 * S1) / np.power(2.0 + n, 4) - (46.617283950617285 * S1) / np.power(2.0 + n, 3) - (42.27160493827161 * S1) / np.power(2.0 + n, 2) - (64.96148967346869 * S1) / (2.0 + n) + (1.7119341563786008 * S2) / (-1.0 + n) - (2.3703703703703702 * S2) / np.power(n, 3) + (0.7407407407407407 * S2) / np.power(n, 2) - (7.703703703703703 * S2) / n + (9.481481481481481 * S2) / np.power(1.0 + n, 3) - (15.012345679012345 * S2) / np.power(1.0 + n, 2) + (27.308641975308642 * S2) / (1.0 + n) - (4.7407407407407405 * S2) / np.power(2.0 + n, 3) + (2.765432098765432 * S2) / np.power(2.0 + n, 2) - (21.020576131687243 * S2) / (2.0 + n) + (0.8559670781893004 * (np.power(S1, 2) + S2)) / (-1.0 + n) - (1.1851851851851851 * (np.power(S1, 2) + S2)) / np.power(n, 3) + (0.37037037037037035 * (np.power(S1, 2) + S2)) / np.power(n, 2) - (2.60082304526749 * (np.power(S1, 2) + S2)) / n + (2.9135802469135803 * (np.power(S1, 2) + S2)) / np.power(1.0 + n, 2) + (5.275720164609053 * (np.power(S1, 2) + S2)) / (1.0 + n) + (2.3703703703703702 * (np.power(S1, 2) + S2)) / np.power(2.0 + n, 3) - (6.518518518518518 * (np.power(S1, 2) + S2)) / np.power(2.0 + n, 2) - (2.8724279835390947 * (np.power(S1, 2) + S2)) / (2.0 + n) - (3.950617283950617 * S21) / n - (2.3703703703703702 * S21) / np.power(1.0 + n, 2) + (7.901234567901234 * S21) / (1.0 + n) + (2.3703703703703702 * S21) / np.power(2.0 + n, 2) - (7.901234567901234 * S21) / (2.0 + n) + (1.1851851851851851 * S211) / n - (2.3703703703703702 * S211) / (1.0 + n) + (2.3703703703703702 * S211) / (2.0 + n) - (3.950617283950617 * S3) / n + (2.3703703703703702 * S3) / np.power(1.0 + n, 2) + (3.1604938271604937 * S3) / (1.0 + n) - (2.3703703703703702 * S3) / np.power(2.0 + n, 2) - (3.1604938271604937 * S3) / (2.0 + n) + (3.950617283950617 * (S1 * S2 - 1.0 * S21 + S3)) / n + (2.3703703703703702 * (S1 * S2 - 1.0 * S21 + S3)) / np.power(1.0 + n, 2) - (7.901234567901234 * (S1 * S2 - 1.0 * S21 + S3)) / (1.0 + n) - (2.3703703703703702 * (S1 * S2 - 1.0 * S21 + S3)) / np.power(2.0 + n, 2) + (7.901234567901234 * (S1 * S2 - 1.0 * S21 + S3)) / (2.0 + n) - (0.6584362139917695 * (np.power(S1, 3) + 3.0 * S1 * S2 + 2.0 * S3)) / n - (0.3950617283950617 * (np.power(S1, 3) + 3.0 * S1 * S2 + 2.0 * S3)) / np.power(1.0 + n, 2) + (1.316872427983539 * (np.power(S1, 3) + 3.0 * S1 * S2 + 2.0 * S3)) / (1.0 + n) + (0.3950617283950617 * (np.power(S1, 3) + 3.0 * S1 * S2 + 2.0 * S3)) / np.power(2.0 + n, 2) - (1.316872427983539 * (np.power(S1, 3) + 3.0 * S1 * S2 + 2.0 * S3)) / (2.0 + n) + (1.1851851851851851 * S31) / n - (2.3703703703703702 * S31) / (1.0 + n) + (2.3703703703703702 * S31) / (2.0 + n) + ( 1.1851851851851851 * ( 0.041666666666666664 * np.power(S1, 4) + 0.25 * np.power(S1, 2) * S2 + 0.125 * np.power(S2, 2) + 0.3333333333333333 * S1 * S3 + 0.25 * S4 ) ) / n - ( 2.3703703703703702 * ( 0.041666666666666664 * np.power(S1, 4) + 0.25 * np.power(S1, 2) * S2 + 0.125 * np.power(S2, 2) + 0.3333333333333333 * S1 * S3 + 0.25 * S4 ) ) / (1.0 + n) + ( 2.3703703703703702 * ( 0.041666666666666664 * np.power(S1, 4) + 0.25 * np.power(S1, 2) * S2 + 0.125 * np.power(S2, 2) + 0.3333333333333333 * S1 * S3 + 0.25 * S4 ) ) / (2.0 + n) + (3.5555555555555554 * S4) / n - (7.111111111111111 * S4) / (1.0 + n) + (7.111111111111111 * S4) / (2.0 + n) - (0.5925925925925926 * (np.power(S2, 2) + S4)) / n + (1.1851851851851851 * (np.power(S2, 2) + S4)) / (1.0 + n) - (1.1851851851851851 * (np.power(S2, 2) + S4)) / (2.0 + n) - (1.1851851851851851 * (S1 * S3 - 1.0 * S31 + S4)) / n + (2.3703703703703702 * (S1 * S3 - 1.0 * S31 + S4)) / (1.0 + n) - (2.3703703703703702 * (S1 * S3 - 1.0 * S31 + S4)) / (2.0 + n) + ( 1.1851851851851851 * (S1 * S21 - 2.0 * S211 + S31 + 0.5 * (np.power(S2, 2) + S4)) ) / n - ( 2.3703703703703702 * (S1 * S21 - 2.0 * S211 + S31 + 0.5 * (np.power(S2, 2) + S4)) ) / (1.0 + n) + ( 2.3703703703703702 * (S1 * S21 - 2.0 * S211 + S31 + 0.5 * (np.power(S2, 2) + S4)) ) / (2.0 + n) - ( 1.1851851851851851 * ( S211 + 0.5 * (S1 * S3 + S1 * (S1 * S2 - 2.0 * S21 + S3) - 2.0 * S31 + S4) ) ) / n + ( 2.3703703703703702 * ( S211 + 0.5 * (S1 * S3 + S1 * (S1 * S2 - 2.0 * S21 + S3) - 2.0 * S31 + S4) ) ) / (1.0 + n) - ( 2.3703703703703702 * ( S211 + 0.5 * (S1 * S3 + S1 * (S1 * S2 - 2.0 * S21 + S3) - 2.0 * S31 + S4) ) ) / (2.0 + n) + (2.5020576131687244 * Sm2) / n + (0.5925925925925926 * Sm2) / np.power(1.0 + n, 2) - (0.6584362139917695 * Sm2) / (1.0 + n) + (1.8436213991769548 * Sm2) / (2.0 + n) - (7.901234567901234 * Sm3) / n + (11.061728395061728 * Sm3) / (1.0 + n) - (11.061728395061728 * Sm3) / (2.0 + n) + (4.7407407407407405 * Sm4) / n - (9.481481481481481 * Sm4) / (1.0 + n) + (9.481481481481481 * Sm4) / (2.0 + n) )
[docs] @nb.njit(cache=True) def gamma_qg_nf1(n, cache, variation): r"""Return the part proportional to :math:`nf^1` of :math:`\gamma_{qg}^{(3)}`. Parameters ---------- n : complex Mellin moment cache: numpy.ndarray Harmonic sum cache variation : int |N3LO| anomalous dimension variation Returns ------- complex |N3LO| non-singlet anomalous dimension :math:`\gamma_{qg}^{(3)}|_{nf^1}` """ S1 = c.get(c.S1, cache, n) S2 = c.get(c.S2, cache, n) S3 = c.get(c.S3, cache, n) S4 = c.get(c.S4, cache, n) S5 = c.get(c.S5, cache, n) common = -7871.5226542038545/np.power(-1. + n,3) + 14103.703703703704/np.power(n,7) + 2588.8395061728397/np.power(n,6) + 68802.34242841466/np.power(n,5) - 35.68779444531073*lm14(n,S1,S2,S3,S4) + 40.511391*lm14m1(n,S1,S2,S3,S4) - 1.8518518518518519*lm15(n,S1,S2,S3,S4,S5) - 2.8806584*lm15m1(n,S1,S2,S3,S4,S5) if variation == 1: fit = -247375.75774004986*(1/(-1. + n) - 1./n) + 56874.46553759208/np.power(-1. + n,2) + 153065.70425758237/np.power(n,4) + 89120.22684116382/np.power(n,3) + 213363.07086259205/np.power(n,2) + 3180.534559480279/n - (33208.154673192555*(S1 - 1.*n*(1.6449340668482262 - 1.*S2)))/np.power(n,2) - 3322.713126819157*lm11(n,S1) + 553.2938548970685*lm12(n,S1,S2) + 158.75214018164516*lm13(n,S1,S2,S3) elif variation == 2: fit = -255300.64168183607*(1/(-1. + n) - 1./n) + 58052.20372466352/np.power(-1. + n,2) + 165947.6629367548/np.power(n,4) + 87940.14555042177/np.power(n,3) + 219563.03597830943/np.power(n,2) + 13782.154548424141/n - (47236.239554654516*(S1 - 1.*n*(1.6449340668482262 - 1.*S2)))/np.power(n,2) + 908.9357647350254*lm12(n,S1,S2) + 171.47745766817872*lm13(n,S1,S2,S3) + 970.4076065981659*lm13m1(n,S1,S2,S3) elif variation == 3: fit = -242660.61714629203*(1/(-1. + n) - 1./n) + 56162.59138894746/np.power(-1. + n,2) + 145568.95877919844/np.power(n,4) + 91620.0763581833/np.power(n,3) + 202918.63844749553/np.power(n,2) + 11377.066521125835/n - (38132.22452502161*(S1 - 1.*n*(1.6449340668482262 - 1.*S2)))/np.power(n,2) + 1026.3816044632167*lm12(n,S1,S2) - 1751.7002864880124*lm12m1(n,S1,S2) + 181.9171385586353*lm13(n,S1,S2,S3) elif variation == 4: fit = -259524.1984624398*(1/(-1. + n) - 1./n) + 58591.28830364899/np.power(-1. + n,2) + 181881.28798856013/np.power(n,4) + 70282.60173332073/np.power(n,3) + 240369.25321380928/np.power(n,2) + 11186.426965436052/n - (81949.00634854026*(S1 - 1.*n*(1.6449340668482262 - 1.*S2)))/np.power(n,2) - 33088.796162944*lm11m1(n,S1) + 1032.4892987335684*lm12(n,S1,S2) + 182.30420495174405*lm13(n,S1,S2,S3) elif variation == 5: fit = -218172.05071206865*(1/(-1. + n) - 1./n) + 52630.53351859396/np.power(-1. + n,2) + 88656.87559909809/np.power(n,4) + 142086.29160312223/np.power(n,3) + 100006.94572043093/np.power(n,2) + 101511.04451818505/n - 90857.77989501068/(1. + n) - (43063.29596432961*(S1 - 1.*n*(1.6449340668482262 - 1.*S2)))/np.power(n,2) + 1068.0436800559826*lm12(n,S1,S2) + 185.87920519454948*lm13(n,S1,S2,S3) elif variation == 6: fit = -228615.46211484703*(1/(-1. + n) - 1./n) + 54086.44766317646/np.power(-1. + n,2) + 122570.71255892045/np.power(n,4) + 91913.78066709182/np.power(n,3) + 198686.11322224894/np.power(n,2) - 21916.281158290127/n - 11188.448464915626*lm11(n,S1) - 288.6043125782358*lm12(n,S1,S2) + 128.6279422497544*lm13(n,S1,S2,S3) - 2297.208338175905*lm13m1(n,S1,S2,S3) elif variation == 7: fit = -279174.74686141574*(1/(-1. + n) - 1./n) + 61675.35655718168/np.power(-1. + n,2) + 203623.91849760094/np.power(n,4) + 72261.1262081702/np.power(n,3) + 283800.62921948434/np.power(n,2) - 52097.21874860658/n - 25731.231535326042*lm11(n,S1) - 2637.2304147864165*lm12(n,S1,S2) + 11813.530413555982*lm12m1(n,S1,S2) + 2.526367648137993*lm13(n,S1,S2,S3) elif variation == 8: fit = -239098.75301911813*(1/(-1. + n) - 1./n) + 55704.755407511606/np.power(-1. + n,2) + 133433.02808202375/np.power(n,4) + 101954.68102915933/np.power(n,3) + 194963.18652126708/np.power(n,2) - 2274.055790556624/n - 5586.551140005561*lm11(n,S1) + 22544.075398836234*lm11m1(n,S1) + 226.80729162572862*lm12(n,S1,S2) + 142.70558679475636*lm13(n,S1,S2,S3) elif variation == 9: fit = -345781.3436260889*(1/(-1. + n) - 1./n) + 71174.93130962232/np.power(-1. + n,2) + 370099.4139304466/np.power(n,4) - 89355.65470903603/np.power(n,3) + 595330.9127886677/np.power(n,2) - 328156.6042167514/n + 306156.82179359894/(1. + n) - 14519.018833525486*lm11(n,S1) - 1181.221350470349*lm12(n,S1,S2) + 67.34403327889338*lm13(n,S1,S2,S3) elif variation == 10: fit = -189717.83410841465*(1/(-1. + n) - 1./n) + 48247.94435904237/np.power(-1. + n,2) + 60212.676453379834/np.power(n,4) + 107033.49639870074/np.power(n,3) + 133203.51842168416/np.power(n,2) + 1303.33515783221/n + 1518.3045318323138*lm12(n,S1,S2) - 9088.7048043846*lm12m1(n,S1,S2) + 225.64383535605478*lm13(n,S1,S2,S3) - 4064.559254773159*lm13m1(n,S1,S2,S3) elif variation == 11: fit = -249553.3220490148*(1/(-1. + n) - 1./n) + 57318.62932737881/np.power(-1. + n,2) + 144265.58335351184/np.power(n,4) + 111968.07357014852/np.power(n,3) + 191250.45508366334/np.power(n,2) + 17314.361281019195/n + 45026.39332999994*lm11m1(n,S1) + 740.8069022955618*lm12(n,S1,S2) + 156.74466224444478*lm13(n,S1,S2,S3) + 2290.9154750749794*lm13m1(n,S1,S2,S3) elif variation == 12: fit = 164981.90805778408*(1/(-1. + n) - 1./n) - 3319.1902098276983/np.power(-1. + n,2) - 708956.7703008109/np.power(n,4) + 700855.4000562581/np.power(n,3) - 1.1337697085772934e6/np.power(n,2) + 1.0068422267096747e6/n - 1.0284779967475284e6/(1. + n) + 2709.9809392914126*lm12(n,S1,S2) + 334.50013143260514*lm13(n,S1,S2,S3) - 10014.262588137737*lm13m1(n,S1,S2,S3) elif variation == 13: fit = -227984.8324412876*(1/(-1. + n) - 1./n) + 54048.98151154018/np.power(-1. + n,2) + 113967.6048415508/np.power(n,4) + 110189.34093808815/np.power(n,3) + 170326.67331718229/np.power(n,2) + 11542.973288722116/n + 28796.031502952967*lm11m1(n,S1) + 1021.0663040321336*lm12(n,S1,S2) - 3276.143876955559*lm12m1(n,S1,S2) + 181.58028948857674*lm13(n,S1,S2,S3) elif variation == 14: fit = -432032.2909414218*(1/(-1. + n) - 1./n) + 83476.22621669284/np.power(-1. + n,2) + 585673.5928852884/np.power(n,4) - 298638.02156309277/np.power(n,3) + 998740.8371967396/np.power(n,2) - 685634.119718282/n + 702608.376556465/(1. + n) + 704.207104145041*lm12(n,S1,S2) - 15297.679777834992*lm12m1(n,S1,S2) + 151.27832196597186*lm13(n,S1,S2,S3) elif variation == 15: fit = -172377.23751488476*(1/(-1. + n) - 1./n) + 46029.384972219035/np.power(-1. + n,2) - 14583.060455699066/np.power(n,4) + 221604.14917968694/np.power(n,3) - 55435.129609225056/np.power(n,2) + 201539.67739043778/n - 191476.85332498729/(1. + n) + 36643.60929648103*lm11m1(n,S1) + 1107.4176978769049*lm12(n,S1,S2) + 189.83827016383634*lm13(n,S1,S2,S3) else: fit = 54050.303305865564/np.power(-1. + n,2) - 228159.14535742637/(-1. + n) + 116361.81262716043/np.power(n,4) + 107389.04759075912/np.power(n,3) + 170221.22878713708/np.power(n,2) + 274171.7610824722/n - 20136.495441164152/(1. + n) - (16239.26140438257*S1)/np.power(n,2) - (16239.26140438257*S2)/n - 4023.197540039458*lm11(n,S1) + 6661.420891021745*lm11m1(n,S1) + 567.3785930765971*lm12(n,S1,S2) - 1173.3798888071458*lm12m1(n,S1,S2) + 164.07463914518564*lm13(n,S1,S2,S3) - 874.313806627577*lm13m1(n,S1,S2,S3) return common + fit
[docs] @nb.njit(cache=True) def gamma_qg_nf2(n, cache, variation): r"""Return the part proportional to :math:`nf^2` of :math:`\gamma_{qg}^{(3)}`. Parameters ---------- n : complex Mellin moment cache: numpy.ndarray Harmonic sum cache variation : int |N3LO| anomalous dimension variation Returns ------- complex |N3LO| non-singlet anomalous dimension :math:`\gamma_{qg}^{(3)}|_{nf^2}` """ S1 = c.get(c.S1, cache, n) S2 = c.get(c.S2, cache, n) S3 = c.get(c.S3, cache, n) S4 = c.get(c.S4, cache, n) S5 = c.get(c.S5, cache, n) common = -1991.111111111111/np.power(n,7) + 2069.3333333333335/np.power(n,6) - 7229.376633440217/np.power(n,5) + 3.511659807956104*lm14(n,S1,S2,S3,S4) - 5.5418381*lm14m1(n,S1,S2,S3,S4) + 0.411522633744856*lm15(n,S1,S2,S3,S4,S5) - 0.8230452700000002*lm15m1(n,S1,S2,S3,S4,S5) if variation == 1: fit = 11302.074334654912*(1/(-1. + n) - 1./n) - 1320.0037968925296/np.power(-1. + n,2) - 11470.302262198706/np.power(n,4) - 1197.1961414701282/np.power(n,3) - 17100.24315576166/np.power(n,2) + 7648.108208393158/n - (1551.5516828229195*(S1 - 1.*n*(1.6449340668482262 - 1.*S2)))/np.power(n,2) + 3557.9066277034826*lm11(n,S1) + 171.03318253716137*lm12(n,S1,S2) - 9.779091112742865*lm13(n,S1,S2,S3) elif variation == 2: fit = 19787.91264185802*(1/(-1. + n) - 1./n) - 2581.106954274633/np.power(-1. + n,2) - 25264.094194264282/np.power(n,4) + 66.41361657930115/np.power(n,3) - 23739.066352155718/np.power(n,2) - 3703.9311863602866/n + (13469.489196307912*(S1 - 1.*n*(1.6449340668482262 - 1.*S2)))/np.power(n,2) - 209.78239089493468*lm12(n,S1,S2) - 23.40515748986745*lm13(n,S1,S2,S3) - 1039.096203888196*lm13m1(n,S1,S2,S3) elif variation == 3: fit = 6253.196339402583*(1/(-1. + n) - 1./n) - 557.7432989957733/np.power(-1. + n,2) - 3442.932697756427/np.power(n,4) - 3873.9952345386782/np.power(n,3) - 5916.536221158395/np.power(n,2) - 1128.6026499620984/n + (3721.064080780992*(S1 - 1.*n*(1.6449340668482262 - 1.*S2)))/np.power(n,2) - 335.54146523765274*lm12(n,S1,S2) + 1875.690344235213*lm12m1(n,S1,S2) - 34.58379657997359*lm13(n,S1,S2,S3) elif variation == 4: fit = 24310.427623197524*(1/(-1. + n) - 1./n) - 3158.3500186622855/np.power(-1. + n,2) - 42325.54228350952/np.power(n,4) + 18973.79622130768/np.power(n,3) - 46018.00103198538/np.power(n,2) - 924.4688406123929/n + (50639.302690081386*(S1 - 1.*n*(1.6449340668482262 - 1.*S2)))/np.power(n,2) + 35430.89974005276*lm11m1(n,S1) - 342.08149811730186*lm12(n,S1,S2) - 34.99826261270785*lm13(n,S1,S2,S3) elif variation == 5: fit = -19968.753127508608*(1/(-1. + n) - 1./n) + 3224.3255608561303/np.power(-1. + n,2) + 57497.57256101853/np.power(n,4) - 57912.364960858366/np.power(n,3) + 104279.57198336071/np.power(n,2) - 97642.54253511738/n + 97288.97428749141/(1. + n) + (9001.170877336448*(S1 - 1.*n*(1.6449340668482262 - 1.*S2)))/np.power(n,2) - 380.1524833607812*lm12(n,S1,S2) - 38.826307723617425*lm13(n,S1,S2,S3) elif variation == 6: fit = 12178.59243809819*(1/(-1. + n) - 1./n) - 1450.2654977988725/np.power(-1. + n,2) - 12895.088918232781/np.power(n,4) - 1066.6754308844847/np.power(n,3) - 17785.97960172896/np.power(n,2) + 6475.534627912801/n + 3190.4037459275073*lm11(n,S1) + 131.69802596523036*lm12(n,S1,S2) - 11.18655302696712*lm13(n,S1,S2,S3) - 107.33021540192573*lm13m1(n,S1,S2,S3) elif variation == 7: fit = 9816.361412874907*(1/(-1. + n) - 1./n) - 1095.6964818770414/np.power(-1. + n,2) - 9108.120591408926/np.power(n,4) - 1984.8871581921073/np.power(n,3) - 13809.258614679915/np.power(n,2) + 5065.420487871078/n + 2510.935618889054*lm11(n,S1) + 21.965477955762765*lm12(n,S1,S2) + 551.9518955699056*lm12m1(n,S1,S2) - 17.078273085286913*lm13(n,S1,S2,S3) elif variation == 8: fit = 11688.792171112142*(1/(-1. + n) - 1./n) - 1374.6549325277456/np.power(-1. + n,2) - 12387.57986068241/np.power(n,4) - 597.544415737091/np.power(n,3) - 17959.922442126444/np.power(n,2) + 7393.258998754254/n + 3452.1356670160476*lm11(n,S1) + 1053.304671083539*lm11m1(n,S1) + 155.7790908183799*lm12(n,S1,S2) - 10.528817402351413*lm13(n,S1,S2,S3) elif variation == 9: fit = 6704.364397603272*(1/(-1. + n) - 1./n) - 651.8568170894251/np.power(-1. + n,2) - 1330.0454328877584/np.power(n,4) - 9535.951065292638/np.power(n,3) + 746.0704525607201/np.power(n,2) - 7832.634543029667/n + 14304.267862325403/(1. + n) + 3034.7927904576964*lm11(n,S1) + 89.99315852559592*lm12(n,S1,S2) - 14.049860818725165*lm13(n,S1,S2,S3) elif variation == 10: fit = 1086.869931803675*(1/(-1. + n) - 1./n) + 214.5933236313328/np.power(-1. + n,2) + 4886.406518103401/np.power(n,4) - 5378.085775495205/np.power(n,3) + 886.4884039541136/np.power(n,2) - 145.57591062121207/n - 383.5448679423884*lm12(n,S1,S2) + 2591.6587045298866*lm12m1(n,S1,S2) - 38.85078659079323*lm13(n,S1,S2,S3) + 396.6326380709743*lm13m1(n,S1,S2,S3) elif variation == 11: fit = 18149.056623377604*(1/(-1. + n) - 1./n) - 2371.9272615951277/np.power(-1. + n,2) - 19081.41595965604/np.power(n,4) - 6785.188814029666/np.power(n,3) - 15665.690406607273/np.power(n,2) - 4711.144806347995/n - 12839.346852387338*lm11m1(n,S1) - 161.84021795598804*lm12(n,S1,S2) - 19.204081107234543*lm13(n,S1,S2,S3) - 1415.6409063305787*lm13m1(n,S1,S2,S3) elif variation == 12: fit = -100056.32530316859*(1/(-1. + n) - 1./n) + 14919.04369686766/np.power(-1. + n,2) + 224216.29758048555/np.power(n,4) - 174707.34699849453/np.power(n,3) + 362165.919896175/np.power(n,2) - 286876.59133906814/n + 293272.1378877304/(1. + n) - 723.3532980336563*lm12(n,S1,S2) - 69.89133322870568*lm13(n,S1,S2,S3) + 2093.200032007698*lm13m1(n,S1,S2,S3) elif variation == 13: fit = 4821.08598119956*(1/(-1. + n) - 1./n) - 351.4904313304405/np.power(-1. + n,2) - 359.1711952380929/np.power(n,4) - 5686.043212666898/np.power(n,3) - 2736.107799606897/np.power(n,2) - 1144.792387449229/n - 2810.008530729486*lm11m1(n,S1) - 335.0227796058678*lm12(n,S1,S2) + 2024.4504536875525*lm12m1(n,S1,S2) - 34.550925607056286*lm13(n,S1,S2,S3) elif variation == 14: fit = 24732.685452531434*(1/(-1. + n) - 1./n) - 3223.0944797587845/np.power(-1. + n,2) - 46389.741636570776/np.power(n,4) + 34208.63219802005/np.power(n,3) - 83575.3975456124/np.power(n,2) + 66887.96728660451/n - 68562.76166626229/(1. + n) - 304.10264672218824*lm12(n,S1,S2) + 3197.550279507657*lm12m1(n,S1,S2) - 31.593962803920494*lm13(n,S1,S2,S3) elif variation == 15: fit = -29540.89833763759*(1/(-1. + n) - 1./n) + 4604.114101824525/np.power(-1. + n,2) + 79077.02907550265/np.power(n,4) - 74533.36724119293/np.power(n,3) + 136770.4440889061/np.power(n,2) - 118550.77103207985/n + 118320.61762162622/(1. + n) - 7659.316849572004*lm11m1(n,S1) - 388.3825020473238*lm12(n,S1,S2) - 39.65383756523093*lm13(n,S1,S2,S3) else: fit = 321.72578082513246/np.power(-1. + n,2) + 84.36283862660154/(-1. + n) + 12108.218046846963/np.power(n,4) - 19333.986960863043/np.power(n,3) + 24036.15277690224/np.power(n,2) - 36952.39876189467/n + 30308.215732860743/(1. + n) + (5018.631677445588*S1)/np.power(n,2) + (5018.631677445588*S2)/n + 1049.7449633329193*lm11(n,S1) + 878.3688118964982*lm11m1(n,S1) - 199.5556809410635*lm12(n,S1,S2) + 682.753445168681*lm12m1(n,S1,S2) - 28.545403117012068*lm13(n,S1,S2,S3) - 4.815643702801845*lm13m1(n,S1,S2,S3) return common + fit
[docs] @nb.njit(cache=True) def gamma_qg(n, nf, cache, variation): r"""Compute the |N3LO| quark-gluon singlet anomalous dimension. Parameters ---------- n : complex Mellin moment nf : int Number of active flavors cache: numpy.ndarray Harmonic sum cache variation : int |N3LO| anomalous dimension variation Returns ------- complex |N3LO| quark-gluon singlet anomalous dimension :math:`\gamma_{qg}^{(3)}(N)` """ return ( +nf * gamma_qg_nf1(n, cache, variation) + nf**2 * gamma_qg_nf2(n, cache, variation) + nf**3 * gamma_qg_nf3(n, cache) )