def pdfflow_tester(pdf, members=None): """Test several pdfflow features: - Check the single/many/all-pid signatures - Checks the python and TF signatures - Checks the output of the python and TF signatures are the same - Checks the expected shape of the signatures is correct """ grid_size = 7 x = np.random.rand(grid_size) q2 = 1000.0 * np.random.rand(grid_size) xtf = float_me(x) q2tf = float_me(q2) # Check I can get just one pid for i in range(-1, 2): # as int res_1 = pdf.py_xfxQ2(i, x, q2) # as list res_2 = pdf.py_xfxQ2([i], x, q2) np.testing.assert_allclose(res_1, res_2) # as tf objects tfpid = int_me([i]) res_3 = pdf.xfxQ2(tfpid, xtf, q2tf) np.testing.assert_allclose(res_2, res_3) # Check shape if members is None: assert res_1.numpy().shape == (grid_size, ) else: assert res_1.numpy().shape == ( members, grid_size, ) # Check I can get more than one pid nfl_size = 6 fl_scheme = pdf.flavor_scheme.numpy() nfl_total = fl_scheme.size many_pid = np.random.choice(fl_scheme, nfl_size) res_1 = pdf.py_xfxQ2(many_pid, x, q2) res_2 = pdf.xfxQ2(int_me(many_pid), xtf, q2tf) np.testing.assert_allclose(res_1, res_2) # Check shape if members is None: assert res_1.numpy().shape == (grid_size, nfl_size) else: assert res_1.numpy().shape == (members, grid_size, nfl_size) # Check I can actually get all PID res_1 = pdf.py_xfxQ2_allpid(x, q2) res_2 = pdf.xfxQ2_allpid(xtf, q2tf) np.testing.assert_allclose(res_1, res_2) # Check shape if members is None: assert res_1.numpy().shape == (grid_size, nfl_total) else: assert res_1.numpy().shape == (members, grid_size, nfl_total)
def py_xfxQ2(self, pid, arr_x, arr_q2): """ Python interface for pdfflow The input gets converted to the right tensorflow type before calling the corresponding function: `xfxQ2` returns PDF evaluated in each f(x,q2) for each pid The output of the function is of shape (members, number_of_points, flavours) but note that dimensions of size 1 will be squeezed out. Example ------- >>> from pdfflow.pflow import mkPDFs >>> from pdfflow.configflow import run_eager, float_me, int_me >>> run_eager(True) >>> pdf = mkPDFs("NNPDF31_nlo_as_0118", [0, 1, 4]) >>> pdf.py_xfxQ2(21, [0.4, 0.5], [15625.0, 15625.0]) <tf.Tensor: shape=(3, 2), dtype=float64, numpy= array([[0.02977381, 0.00854525], [0.03653673, 0.00929325], [0.031387 , 0.00896622]])> >>> pdf.py_xfxQ2([1,2], [0.4], [15625.0]) <tf.Tensor: shape=(3, 2), dtype=float64, numpy= array([[0.05569674, 0.19323399], [0.05352555, 0.18965438], [0.04515956, 0.18704451]])> Parameters ---------- pid: list(int) list of PID to be computed arr_x: np.array grid on x where to compute the PDF arr_q2: np.array grid on q^2 where to compute the PDF Returns ------- pdf: tensor grid of results ([members], [number of points], [flavour]) """ # if pid is None, the mask is set to true everywhere if pid is None: return self.py_xfxQ2_allpid(arr_x, arr_q2) tensor_pid = tf.reshape(int_me(pid), (-1, )) a_x = float_me(arr_x) a_q2 = float_me(arr_q2) return self.xfxQ2(tensor_pid, a_x, a_q2)
def py_alphasQ(self, arr_q): """ User interface for pdfflow alphas interpolation when called with python variables. Returns alpha_s for each value of Q User interface for pdfflow alphas interpolation when called with Example ------- >>> from pdfflow.pflow import mkPDF >>> from pdfflow.configflow import run_eager >>> run_eager(True) >>> pdf = mkPDF("NNPDF31_nlo_as_0118/0") >>> pdf.py_alphasQ([125.0, 94]) <tf.Tensor: shape=(2,), dtype=float64, numpy=array([0.11264939, 0.11746421])> Parameters ---------- arr_q: tf.tensor, dtype=float grid on q where to compute alphas Returns ------- alphas: tensor alphas evaluated in each q query point """ # Parse the input a_q = float_me(arr_q) # Perform the actual computation return self.alphasQ(a_q)
def rapidity_dif(p1, p2): """ Computes the rapidity difference between p1 and p2 y = 1/2*log(E+z / E-z) """ num = (p1[0, :] + p1[3, :]) * (p2[0, :] - p2[3, :]) den = (p1[0, :] - p1[3, :]) * (p2[0, :] + p2[3, :]) return float_me(0.5) * tf.math.log(num / den)
def lowx_lowq2_extrapolation( a_x, a_q2, padded_x, s_x, padded_q2, s_q2, actual_padded, ): """ Extrapolation in low q2, low x regime Parameters ---------- a_x: tf.tensor of shape [None] query of values of log(x) a_q2: tf.tensor of shape [None] query of values of log(q2) padded_x: tf.tensor of shape [None] value for all the knots on the x axis padded with one zero at the beginning and one at the end to avoid out of range errors when queryingpoints near boundaries s_x: tf.tensor of shape [] size of x knots tensor without padding padded_q2: tf.tensor of shape [None] value for all the knots on the q2 axis padded with one zero at the beginning and one at the end to avoid out of range errors when querying points near boundaries s_q2: tf.tensor of shape [] size of q2 knots tensor without padding actual_padded: tf.tensor of shape [None,None] pdf values: first axis is the flattened padded (q2,x) grid, second axis is needed pid column (dimension depends on the query) """ corn_x = padded_x[1:3] corn_q2 = tf.stack([padded_q2[1], padded_q2[1]], 0) f = interpolate( tf.concat([corn_x, corn_x], 0), tf.concat([corn_q2, 1.01 * corn_q2], 0), padded_x, s_x, padded_q2, s_q2, actual_padded, ) fq2Min = extrapolate_linear(a_x, corn_x[0], corn_x[1], f[:1], f[1:2]) fq2Min1 = extrapolate_linear(a_x, corn_x[0], corn_x[1], f[2:3], f[3:]) a_q2 = tf.expand_dims(tf.math.exp(a_q2), 1) corn_q2 = tf.math.exp(corn_q2[0]) mask = tf.math.abs(fq2Min) >= 1e-5 anom = tf.where( mask, tf.maximum(float_me(-2.5), (fq2Min1 - fq2Min) / fq2Min / 0.01), fone ) factor = tf.math.pow(a_q2 / corn_q2, anom * a_q2 / corn_q2 + 1.0 - a_q2 / corn_q2) return fq2Min * factor
def py_xfxQ2_allpid(self, arr_x, arr_q2): """ Python interface for pdfflow The input gets converted to the right tensorflow type before calling the corresponding function: `xfxQ2_allpid` Returns all flavours The output of the function is of shape (members, number_of_points, all_flavours) but note that dimensions of size 1 will be squeezed out. Example ------- >>> from pdfflow.pflow import mkPDF >>> from pdfflow.configflow import run_eager >>> run_eager(True) >>> pdf = mkPDF("NNPDF31_nlo_as_0118/0") >>> pdf.py_xfxQ2_allpid([0.4], [15625.0]) <tf.Tensor: shape=(11,), dtype=float64, numpy= array([2.03487396e-04, 4.64913697e-03, 2.36497526e-04, 4.54391659e-03, 3.81215383e-03, 6.68418154e-02, 8.91427455e-02, 2.86122760e-01, 5.96581806e-03, 4.64913709e-03, 2.03487286e-04])> Parameters ---------- arr_x: np.array grid on x where to compute the PDF arr_q2: np.array grid on q^2 where to compute the PDF Returns ------- pdf: tensor grid of results ([members], [number of points], flavour) """ a_x = float_me(arr_x) a_q2 = float_me(arr_q2) return self.xfxQ2_allpid(a_x, a_q2)
prop = propagator_w(sa1) * propagator_w(sb2) coup = tf.square(mw / tf.pow(stw, 1.5)) rmcom = coup / prop # Compute the amplitude # W-boson, so only Left-Left amplitude = zA(pa, pb) * zB(p1, p2) amp2 = tf.math.real(amplitude * tf.math.conj(amplitude)) me_res = 2.0 * amp2 * rmcom return me_res # Leading Order matrix element factor_lo = float_me(1.0702411577062499e-4) # there is no alpha_s, alpha_ew computed at Mz val @tf.function(input_signature=4 * [TFLOAT4]) def qq_h_lo(pa, pb, p1, p2): """ Computes the LO q Q -> Q q H (WW->H) """ return factor_lo * partial_lo(pa, pb, p1, p2) @tf.function(input_signature=[TFLOAT4] * 5) def partial_qq_h_qQg(pa, pb, p1, p2, p3): """ Gluon radiated from leg pa-p1 """ pa13 = pa - (p1 + p3) sa13 = dot_product(pa13, pa13) sb2 = -2.0 * dot_product(pb, p2) prop = propagator_w(sa13) * propagator_w(sb2)
def __init__(self, grid, i=0, total=0, compile_functions=True, alpha_s=False): name_sg = f"grid_{i}" self.alpha_s = alpha_s if alpha_s: name_sg += "_alpha" self.name_sg = name_sg super().__init__(name=f"Parent_{name_sg}") q2min = min(grid.q2) q2max = max(grid.q2) self.log_q2min = float_me(np.log(q2min)) self.log_q2max = float_me(np.log(q2max)) # Save grid shape information self.s_q2 = int_me(grid.q2.size) # Insert a padding at the beginning and the end log_q2pad = np.pad(np.log(grid.q2), 1, mode="edge") log_q2pad[0] *= 0.99 log_q2pad[-1] *= 1.01 self.padded_q2 = float_me(log_q2pad) # Depending on whether it is an alphs_s grid or a pdf grid # we might need to change some options compilation_options = OPT.copy() if alpha_s: # the grid is sized (q.size), pad it with 0s self.padded_grid = float_me(np.pad(grid.grid, (1, 1))) if i == 0: self.fn_interpolation = alphas_first_subgrid elif i == (total - 1): self.fn_interpolation = alphas_last_subgrid else: self.fn_interpolation = alphas_inner_subgrid # Change the function signature to that of alpha_s compilation_options[ "input_signature"] = ALPHAS_GRID_FUNCTION_SIGNATURE else: # If this is a pdf grid, save also the x information xmin = min(grid.x) xmax = max(grid.x) self.log_xmin = float_me(np.log(xmin)) self.log_xmax = float_me(np.log(xmax)) self.s_x = int_me(grid.x.size) log_xpad = np.pad(np.log(grid.x), 1, mode="edge") log_xpad[0] *= 0.99 log_xpad[-1] *= 1.01 self.padded_x = float_me(log_xpad) # Finally parse the grid # the grid is sized (x.size * q.size, flavours) reshaped_grid = grid.grid.reshape(grid.x.size, grid.q2.size, -1) # and pad it with 0s in x and q padded_grid = np.pad(reshaped_grid, ((1, 1), (1, 1), (0, 0))) # flatten the x and q dimensions again and store it self.padded_grid = float_me(padded_grid.reshape( -1, grid.flav.size)) # Depending on the index of the grid, select with interpolation function should be run if i == 0: self.fn_interpolation = first_subgrid elif i == (total - 1): self.fn_interpolation = last_subgrid else: self.fn_interpolation = inner_subgrid if compile_functions: self.fn_interpolation = tf.function(self.fn_interpolation, **compilation_options)
def main(pdfname, pid, no_tex=True): """Testing PDFflow vs LHAPDF performance.""" mpl.rcParams['text.usetex'] = no_tex mpl.rcParams['savefig.format'] = 'pdf' mpl.rcParams['figure.figsize'] = [11, 5.5] mpl.rcParams['axes.titlesize'] = 20 mpl.rcParams['ytick.labelsize'] = 17 mpl.rcParams['xtick.labelsize'] = 17 mpl.rcParams['legend.fontsize'] = 14 import pdfflow.pflow as pdf p = pdf.mkPDF(pdfname, DIRNAME) l_pdf = lhapdf.mkPDF(pdfname) name = '\_'.join(pdfname.split('_')) s = time.time() p.trace() print("\nPDFflow\n\tBuilding graph time: %f\n" % (time.time() - s)) fig = plt.figure() gs = fig.add_gridspec(nrows=1, ncols=2, wspace=0.05) ax = fig.add_subplot(gs[0]) x = np.logspace(-12, 0, 10000, dtype=float) q2 = np.array([1.65, 1.7, 4.92, 1e2, 1e3, 1e4, 1e5, 1e6, 2e6], dtype=float)**2 for iq2 in q2: vl = np.array([l_pdf.xfxQ2(pid, ix, iq2) for ix in x]) vp = p.py_xfxQ2(pid, float_me(x), float_me([iq2] * len(x))) ax.plot(x, np.abs(vp - vl) / (np.abs(vl) + EPS), label=r'$Q=%s$' % sci_notation(iq2**0.5, 2)) ax.hlines(1e-3, plt.xlim()[0], plt.xlim()[1], linestyles='dotted', color='red') ax.set_xscale('log') ax.set_yscale('log') ax.set_xlim([1e-12, 1.]) ax.set_ylim([EPS, .01]) ax = set_ticks(ax, -12, 0, 13, 'x', 4) ax.tick_params(axis='x', which='both', direction='in', bottom=True, labelbottom=True, top=True, labeltop=False) ax = set_ticks(ax, -15, -3, 16, 'y') ax.tick_params(axis='y', which='both', direction='in', left=True, labelleft=True, right=True, labelright=False) ax.set_title(r'%s, flav = %d' % (name, pid)) ylabel = r'$\displaystyle{r_{i}(x,Q)}$' if no_tex else '$r_{i}$(x,Q)' ax.set_ylabel(ylabel, fontsize=20) ax.set_xlabel(r'$x$', fontsize=17) ax.legend(frameon=False, ncol=2, loc='upper right', bbox_to_anchor=(1.02, 0.9)) x = np.array([1e-10, 1e-9, 1.1e-9, 5e-7, 1e-6, 1e-4, 1e-2, 0.5, 0.99], dtype=float) q2 = np.logspace(1, 7, 10000, dtype=float)**2 ax = fig.add_subplot(gs[1]) for ix in x: s_time = time.time() vl = np.array([l_pdf.xfxQ2(pid, ix, iq2) for iq2 in q2]) l_time = time.time() vp = p.py_xfxQ2(pid, float_me([ix] * len(q2)), float_me(q2)) p_time = time.time() ax.plot(q2**0.5, np.abs(vp - vl) / (np.abs(vl) + EPS), label=r'$x=%s$' % sci_notation(ix, 1)) ax.hlines(1e-3, plt.xlim()[0], plt.xlim()[1], linestyles='dotted', color='red') ax.set_xscale('log') ax.set_yscale('log') ax.set_xlim([1, 1e7]) ax.set_ylim([EPS, .01]) ax = set_ticks(ax, 1, 7, 9, 'x') ax.tick_params(axis='x', which='both', direction='in', top=True, labeltop=False, bottom=True, labelbottom=True) ax = set_ticks(ax, -15, -3, 16, 'y') ax.tick_params(axis='y', which='both', direction='in', right=True, labelright=False, left=True, labelleft=False) ax.set_title(r'%s, flav = %d' % (name, pid)) ax.set_xlabel(r'$Q$', fontsize=17) ax.legend(frameon=False, ncol=2, loc='upper right', bbox_to_anchor=(1.02, 0.9)) plt.savefig('diff_%s_flav%d.pdf' % (pdfname.replace('/', '-'), pid), bbox_inches='tight', dpi=250) plt.close() print("\nDry run time comparison:") print("{:>10}:{:>15.8f}".format("lhapdf", l_time - s_time)) print("{:>10}:{:>15.8f}".format("pdfflow", p_time - l_time))
def lowq2_extrapolation( a_x, a_q2, padded_x, s_x, padded_q2, s_q2, actual_padded, ): """ Extrapolation in low q2 regime Parameters ---------- a_x: tf.tensor of shape [None] query of values of log(x) a_q2: tf.tensor of shape [None] query of values of log(q2) padded_x: tf.tensor of shape [None] value for all the knots on the x axis padded with one zero at the beginning and one at the end to avoid out of range errors when queryingpoints near boundaries s_x: tf.tensor of shape [] size of x knots tensor without padding padded_q2: tf.tensor of shape [None] value for all the knots on the q2 axis padded with one zero at the beginning and one at the end to avoid out of range errors when querying points near boundaries s_q2: tf.tensor of shape [] size of q2 knots tensor without padding actual_padded: tf.tensor of shape [None,None] pdf values: first axis is the flattened padded (q2,x) grid, second axis is needed pid column (dimension depends on the query) """ corn_q2 = tf.stack([padded_q2[1], 1.01 * padded_q2[1]], 0) x, q2 = tf.meshgrid(a_x, corn_q2) s = tf.size(a_x, out_type=DTYPEINT) fq2Min = interpolate( tf.reshape(x, [-1]), tf.reshape(q2, [-1]), padded_x, s_x, padded_q2, s_q2, actual_padded, ) fq2Min1 = fq2Min[s:] fq2Min = fq2Min[:s] a_q2 = tf.math.exp(a_q2) corn_q2 = tf.math.exp(corn_q2[:1]) mask = tf.math.abs(fq2Min) >= 1e-5 anom = tf.where( mask, tf.maximum(float_me(-2.5), (fq2Min1 - fq2Min) / fq2Min / 0.01), fone ) corn_q2 = tf.expand_dims(corn_q2, 1) a_q2 = tf.expand_dims(a_q2, 1) return fq2Min * tf.math.pow( a_q2 / corn_q2, anom * a_q2 / corn_q2 + 1.0 - a_q2 / corn_q2 )
def test_float_me(): res = float_me(4.0) assert res.dtype == DTYPE
TECH_CUT, shat_min, s_in, higgs_mass, pt2_cut, rdistance, deltaycut, m2jj_cut, ) import spinors # Control flags UNIT_PHASE = False # Constants tfpi = float_me(3.1415926535897932385) costhmax = fone costhmin = float_me(-1.0) * fone phimin = fzero phimax = float_me(2.0) * tfpi # Jet separation @tf.function def rapidity_dif(p1, p2): """ Computes the rapidity difference between p1 and p2 y = 1/2*log(E+z / E-z) """ num = (p1[0, :] + p1[3, :]) * (p2[0, :] - p2[3, :]) den = (p1[0, :] - p1[3, :]) * (p2[0, :] + p2[3, :]) return float_me(0.5) * tf.math.log(num / den)