Skip to content

A set of python-based Monte Carlo tools for jet physics, including Monte Carlo integration, parton showering, and examples.

License

Notifications You must be signed in to change notification settings

samcaf/JetMonteCarlo

Repository files navigation

JetMC

standard-readme compliant

A set of python-based Monte Carlo tools for jet physics, including Monte Carlo integration, parton showering, and examples.

The tools in this repository are designed in particular to produce distributions of correlations in jet substructure at leading order (LO), leading logarithmic (LL), and modified leading logarithmic (MLL) accuracy. This repository contains tools which allow two separate, monte carlo based approaches to the calculation of substructure correlations: Monte Carlo integration over the phase space of particles in jets, and parton showering algorithms. Each of these are described in more detail below.

Table of Contents

Monte Carlo Integration

Monte Carlo integration is a method of sampling over a parameter space in order to numerically estimate the integrals of distributions over that parameter space. This repository includes some very basic tools for general Monte Carlo integration which will guide us toward our goal of jet substructure. We also include some simple examples which use these tools.

The more specific goal of this repository, however, is to sample over the phase space of emissions within jets. The jet samplers in this repository sample over the leading logarithmic phase space of ungroomed jets, as well as jets which are groomed through several complicated procedures. This phase space is weighted by the probability distributions and psuedo-probability distributions which are derived from perturbative Quantum Chromodynamics (QCD). In order to characterize the substructure of jets, we are interested in using Monte Carlo integration to obtain probability distribution functions (pdfs) and cumulative distribution functions (cdfs) of correlations between particles within a jet.

Incomplete

  • Integrate 2D
    • Integrate 2D Tests

Parton Showering Techniques

Parton showering is another way to perform calculations which will help us characterize jet substructure. The procedure of parton showering is, in broad strokes, a method of producing a tree of parton splittings which emulate the structure of perturbative QCD. These parton splittings describe the splitting of a single parton into two daughter partons, and they will be ordered by the value of some observable (i.e. the splittings are ordered by some analogue of "time" in the parton shower). After the splittings reach a scale which we expect will be dominated by non-perturbative dynamics, rather than perturbative QCD, we stop splitting the partons.

Our parton showers are ordered by the observable angularity. The angularity of each subsequent splitting in our parton shower is generated by using the inverse transform method using leading logarithmic distributions for angularity. The splittings stop when the angularity reaches a non-perturbative scale of 1 GeV. We can also include effects at higher orders of accuracy, such as effects due to non-singular pieces of splitting functions and effects due to the running coupling of QCD. We do this by modifying the above procedure with an additional stage of acceptance-rejection sampling, which allows us to reproduce the distributions including these non-trivial effects by starting from leading logarithmic accuracy.

Usage

Basic Monte Carlo Integration

The fundamental tools of Monte Carlo integration are random number generators which sample over a phase space of interest (samplers). A simple example is a sampler which samples uniformly over the unit interval (0,1):

# ------------------------------------
# Simple sampler class:
# ------------------------------------
class simpleSampler(sampler):
    """A very simple linspace sampler, which samples from 0 to 1."""
    # ------------------
    # Sampling:
    # ------------------
    def setArea(self):
        if self.sampleMethod=='lin':
            self.area = 1.
        elif self.sampleMethod=='log':
            self.area = np.log(1./self.epsilon)

    def addSample(self):
        if self.sampleMethod=='lin':
            sample = getLinSample(0., 1.)
            jac    = 1.
        if self.sampleMethod=='log':
            sample = getLogSample(0., 1., self.epsilon)
            jac    = sample
        self.samples.append(sample)
        self.jacobians.append(jac)

    # ------------------
    # Init:
    # ------------------
    def __init__(self, sampleMethod, **kwargs):
        super().__init__(sampleMethod,**kwargs)

After sampling, we integrate a function over the desired phase space to obtain a numerical integral. We introduce integrators to do this for us. The integrator class integrates over the phase space to obtain the integral of a given function over subregions of the parameter space, as well as the whole parameter space. As an example, we can integrate simple polynomials over the phase space of our simple sampler

# Sampling
testSampler = simpleSampler('lin')
testSampler.generateSamples(numSamples)
samples = testSampler.getSamples()

# Setting up integrator
testInt = integrator()
testInt.setFirstBinBndCondition(0.)
# Telling the integrator that the initial value of the integral should be zero.
testInt.setBins(numBins, samples, 'lin')
# Giving the integrator several bins.
# The upper/right limits of these bins will be the upper bounds of our integrations

# Obtaining weights (corresponding to the integrand), binned observables, and area of the integration region
weights = samples**n / (n+1)
jacs    = testSampler.jacobians
obs     = samples
area    = testSampler.area

# Setting the integrand, and obtaining the integral
testInt.setDensity(obs, weights * jacs, area)
testInt.integrate()

# Getting the value of the integral, the associated area, and associated x values
integral = testInt.integral
intErr   = testInt.integralErr
xs       = testInt.bins[1:]

This gives us precisely the integrals we expect:

simple_int

Basic Jet Monte Carlo Integration

When performing integration over the phase space of emissions within a jet, the procedure is almost identical. The only complication is that this phase space is dominated by infra-red, low energy emissions. When formulated precisely, this phenomenon is described through infra-red and collinear singularities in the probability distributions for quark and gluon splittings.

Some simple examples are included in our tests folder. Our procedures are able to numerically reproduce the analytic expressions for more complicated observables. Here is an example in which we use Monte Carlo integration to reproduce the probability distribution for the energy and angle of emissions within a groomed gluon jet:

LL_pdf

Basic Parton Showering

We can implement a basic parton shower very easily:

from jetmontecarlo.utils.partonshower_utils import *

# Setting up initial parameters for the jet: radius, initial angularity, and initial transverse momentum:
radius = 1
ang_init = radius**beta / 2.
momentum = Vector([0,PT,0])

# Setting up the jet
jet      = Jet(momentum, radius, partons=None)
mother   = jet.partons[0]

# Performing the parton shower
angularity_shower(mother, ang_init, beta, 'quark', jet,
                      acc=acc, split_soft=False)
jet.has_showered = True

The soul of the parton shower is in the angularity_shower method. This can be roughly described with the following code flowchart:

ps_codeflow

Here is a characterisctic visualization of the corresponding tree of showered partons in our jet, with a green cone to indicate the jet cone and green lines and dots to indicate final state partons with no further splittings:

shower_vis

After producing a set of parton showers representing jets at leading logarithmic or higher levels of accuracy, we are ready to calculate the distributions of jet substructure observables. Our simplest example is a bit too long to include here, but produces results which agree with analytic leading logarithmic distributions. Here is an example showing that our parton shower accurately reproduces quark jet angularities at leading logarithmic accuracy:

LL_partonshower

To Add

Additional Utilities

I would like to add additional utilities, including observables, formatting utilities (using common formatting for events such as in .hepmc files), and bridges to other technology (such as FastJet).

Multivariable Parton Showering

I'd like to add functionality for parton showers with multiple evolution variables, as in

https://arxiv.org/pdf/2201.08056.pdf

Evolution Equations

I'd like to add additonal utilities for solving evolution equations, both analytically and numerically. Basic equations include DGLAP equations, and perhaps BFKL equations. Furthermore, it seems like parton branching is a place where using Monte Carlo to solve evolution equations has been well explored, and provides transverse momentum dependent (TMD) information:

https://arxiv.org/pdf/1708.03279.pdf

https://arxiv.org/pdf/1704.01757.pdf

https://arxiv.org/pdf/2202.00969.pdf

Maintainers

@samcaf

Contributing

PRs accepted.

Small note: If editing the README, please conform to the standard-readme specification.

License

MIT © 2021 Samuel Alipour-fard

About

A set of python-based Monte Carlo tools for jet physics, including Monte Carlo integration, parton showering, and examples.

Topics

Resources

License

Stars

Watchers

Forks

Releases

No releases published

Packages

No packages published