madminer.delphes module

class madminer.delphes.DelphesReader(filename)[source]

Bases: object

Detector simulation with Delphes and simple calculation of observables.

After setting up the parameter space and benchmarks and running MadGraph and Pythia, all of which is organized in the madminer.core.MadMiner class, the next steps are the simulation of detector effects and the calculation of observables. Different tools can be used for these tasks, please feel free to implement the detector simulation and analysis routine of your choice.

This class provides an example implementation based on Delphes. Its workflow consists of the following steps:

  • Initializing the class with the filename of a MadMiner HDF5 file (the output of madminer.core.MadMiner.save())
  • Adding one or multiple event samples produced by MadGraph and Pythia in DelphesProcessor.add_sample().
  • Running Delphes on the samples that require it through DelphesProcessor.run_delphes().
  • Optionally, acceptance cuts for all visible particles can be defined with DelphesProcessor.set_acceptance().
  • Defining observables through DelphesProcessor.add_observable() or DelphesProcessor.add_observable_from_function(). A simple set of default observables is provided in DelphesProcessor.add_default_observables()
  • Optionally, cuts can be set with DelphesProcessor.add_cut()
  • Calculating the observables from the Delphes ROOT files with DelphesProcessor.analyse_delphes_samples()
  • Saving the results with DelphesProcessor.save()

Please see the tutorial for a detailed walk-through.

Parameters:
filename : str or None, optional

Path to MadMiner file (the output of madminer.core.MadMiner.save()). Default value: None.

Methods

add_cut(self, definition[, pass_if_not_parsed]) Adds a cut as a string that can be parsed by Python’s eval() function and returns a bool.
add_default_observables(self[, …]) Adds a set of simple standard observables: the four-momenta (parameterized as E, pT, eta, phi) of the hardest visible particles, and the missing transverse energy.
add_observable(self, name, definition[, …]) Adds an observable as a string that can be parsed by Python’s eval() function.
add_observable_from_function(self, name, fn) Adds an observable defined through a function.
add_sample(self, hepmc_filename, …[, …]) Adds a sample of simulated events.
analyse_delphes_samples(self[, …]) Main function that parses the Delphes samples (ROOT files), checks acceptance and cuts, and extracts the observables and weights.
reset_cuts(self) Resets all cuts.
reset_observables(self) Resets all observables.
run_delphes(self, delphes_directory, …[, …]) Runs the fast detector simulation Delphes on all HepMC samples added so far for which it hasn’t been run yet.
save(self, filename_out[, shuffle]) Saves the observable definitions, observable values, and event weights in a MadMiner file.
set_acceptance(self[, pt_min_e, pt_min_mu, …]) Sets acceptance cuts for all visible particles.
add_cut(self, definition, pass_if_not_parsed=False)[source]

Adds a cut as a string that can be parsed by Python’s eval() function and returns a bool.

Parameters:
definition : str

An expression that can be parsed by Python’s eval() function and returns a bool: True for the event to pass this cut, False for it to be rejected. In the definition, all visible particles can be used: e, mu, j, a, and l provide lists of electrons, muons, jets, photons, and leptons (electrons and muons combined), in each case sorted by descending transverse momentum. met provides a missing ET object. visible and all provide access to the sum of all visible particles and the sum of all visible particles plus MET, respectively. All these objects are instances of MadMinerParticle, which inherits from scikit-hep’s [LorentzVector](http://scikit-hep.org/api/math.html#vector-classes). See the link for a documentation of their properties. In addition, MadMinerParticle have properties charge and pdg_id, which return the charge in units of elementary charges (i.e. an electron has e[0].charge = -1.), and the PDG particle ID. For instance, “len(e) >= 2” requires at least two electrons passing the acceptance cuts, while “mu[0].charge > 0.” specifies that the hardest muon is positively charged.

pass_if_not_parsed : bool, optional

Whether the cut is passed if the observable cannot be parsed. Default value: False.

Returns:
None
add_default_observables(self, n_leptons_max=2, n_photons_max=2, n_jets_max=2, include_met=True, include_visible_sum=True, include_numbers=True, include_charge=True)[source]

Adds a set of simple standard observables: the four-momenta (parameterized as E, pT, eta, phi) of the hardest visible particles, and the missing transverse energy.

Parameters:
n_leptons_max : int, optional

Number of hardest leptons for which the four-momenta are saved. Default value: 2.

n_photons_max : int, optional

Number of hardest photons for which the four-momenta are saved. Default value: 2.

n_jets_max : int, optional

Number of hardest jets for which the four-momenta are saved. Default value: 2.

include_met : bool, optional

Whether the missing energy observables are stored. Default value: True.

include_visible_sum : bool, optional

Whether observables characterizing the sum of all particles are stored. Default value: True.

include_numbers : bool, optional

Whether the number of leptons, photons, and jets is saved as observable. Default value: True.

include_charge : bool, optional

Whether the lepton charge is saved as observable. Default value: True.

Returns:
None
add_observable(self, name, definition, required=False, default=None)[source]

Adds an observable as a string that can be parsed by Python’s eval() function.

Parameters:
name : str

Name of the observable. Since this name will be used in eval() calls for cuts, this should not contain spaces or special characters.

definition : str

An expression that can be parsed by Python’s eval() function. As objects, the visible particles can be used: e, mu, j, a, and l provide lists of electrons, muons, jets, photons, and leptons (electrons and muons combined), in each case sorted by descending transverse momentum. met provides a missing ET object. visible and all provide access to the sum of all visible particles and the sum of all visible particles plus MET, respectively. All these objects are instances of MadMinerParticle, which inherits from scikit-hep’s [LorentzVector](http://scikit-hep.org/api/math.html#vector-classes). See the link for a documentation of their properties. In addition, MadMinerParticle have properties charge and pdg_id, which return the charge in units of elementary charges (i.e. an electron has e[0].charge = -1.), and the PDG particle ID. For instance, “abs(j[0].phi() - j[1].phi())” defines the azimuthal angle between the two hardest jets.

required : bool, optional

Whether the observable is required. If True, an event will only be retained if this observable is successfully parsed. For instance, any observable involving “j[1]” will only be parsed if there are at least two jets passing the acceptance cuts. Default value: False.

default : float or None, optional

If required=False, this is the placeholder value for observables that cannot be parsed. None is replaced with np.nan. Default value: None.

Returns:
None
add_observable_from_function(self, name, fn, required=False, default=None)[source]

Adds an observable defined through a function.

Parameters:
name : str

Name of the observable. Since this name will be used in eval() calls for cuts, this should not contain spaces or special characters.

fn : function

A function with signature observable(leptons, photons, jets, met) where the input arguments are lists of MadMinerParticle instances and a float is returned. The function should raise a RuntimeError to signal that it is not defined.

required : bool, optional

Whether the observable is required. If True, an event will only be retained if this observable is successfully parsed. For instance, any observable involving “j[1]” will only be parsed if there are at least two jets passing the acceptance cuts. Default value: False.

default : float or None, optional

If required=False, this is the placeholder value for observables that cannot be parsed. None is replaced with np.nan. Default value: None.

Returns:
None
add_sample(self, hepmc_filename, sampled_from_benchmark, is_background=False, delphes_filename=None, lhe_filename=None, k_factor=1.0, weights='lhe', systematics=None)[source]

Adds a sample of simulated events. A HepMC file (from Pythia) has to be provided always, since some relevant information is only stored in this file. The user can optionally provide a Delphes file, in this case run_delphes() does not have to be called.

By default, the weights are read out from the Delphes file and their names from the HepMC file. There are some issues with current MadGraph versions that lead to Pythia not storing the weights. As work-around, MadMiner supports reading weights from the LHE file (the observables still come from the Delphes file). To enable this, use weights=”lhe”.

Parameters:
hepmc_filename : str

Path to the HepMC event file (with extension ‘.hepmc’ or ‘.hepmc.gz’).

sampled_from_benchmark : str

Name of the benchmark that was used for sampling in this event file (the keyword sample_benchmark of madminer.core.MadMiner.run()).

is_background : bool, optional

Whether the sample is a background sample (i.e. without benchmark reweighting).

delphes_filename : str or None, optional

Path to the Delphes event file (with extension ‘.root’). If None, the user has to call run_delphes(), which will create this file. Default value: None.

lhe_filename : None or str, optional

Path to the LHE event file (with extension ‘.lhe’ or ‘.lhe.gz’). This is only needed if weights is “lhe”.

k_factor : float, optional

Multiplies the cross sections found in the sample. Default value: 1.

weights : {“delphes”, “lhe”}, optional

If “delphes”, the weights are read out from the Delphes ROOT file, and their names are taken from the HepMC file. If “lhe” (and lhe_filename is not None), the weights are taken from the LHE file (and matched with the observables from the Delphes ROOT file). The “delphes” behaviour is generally better as it minimizes the risk of mismatching observables and weights, but for some MadGraph and Delphes versions there are issues with weights not being saved in the HepMC and Delphes ROOT files. In this case, setting weights to “lhe” and providing the unweighted LHE file from MadGraph may be an easy fix. Default value: “lhe”.

systematics : None or list of str, optional

List of systematics associated with this sample. Default value: None.

Returns:
None
analyse_delphes_samples(self, generator_truth=False, delete_delphes_files=False, reference_benchmark=None, parse_lhe_events_as_xml=True)[source]

Main function that parses the Delphes samples (ROOT files), checks acceptance and cuts, and extracts the observables and weights.

Parameters:
generator_truth : bool, optional

If True, the generator truth information (as given out by Pythia) will be parsed. Detector resolution or efficiency effects will not be taken into account.

delete_delphes_files : bool, optional

If True, the Delphes ROOT files will be deleted after extracting the information from them. Default value: False.

reference_benchmark : str or None, optional

The weights at the nuisance benchmarks will be rescaled to some reference theta benchmark: dsigma(x|theta_sampling(x),nu) -> dsigma(x|theta_ref,nu) = dsigma(x|theta_sampling(x),nu) * dsigma(x|theta_ref,0) / dsigma(x|theta_sampling(x),0). This sets the name of the reference benchmark. If None, the first one will be used. Default value: None.

parse_lhe_events_as_xml : bool, optional

Decides whether the LHE events are parsed with an XML parser (more robust, but slower) or a text parser (less robust, faster). Default value: True.

Returns:
None
reset_cuts(self)[source]

Resets all cuts.

reset_observables(self)[source]

Resets all observables.

run_delphes(self, delphes_directory, delphes_card, initial_command=None, log_file=None)[source]

Runs the fast detector simulation Delphes on all HepMC samples added so far for which it hasn’t been run yet.

Parameters:
delphes_directory : str

Path to the Delphes directory.

delphes_card : str

Path to a Delphes card.

initial_command : str or None, optional

Initial bash commands that have to be executed before Delphes is run (e.g. to load the correct virtual environment). Default value: None.

log_file : str or None, optional

Path to log file in which the Delphes output is saved. Default value: None.

Returns:
None
save(self, filename_out, shuffle=True)[source]

Saves the observable definitions, observable values, and event weights in a MadMiner file. The parameter, benchmark, and morphing setup is copied from the file provided during initialization. Nuisance benchmarks found in the HepMC file are added.

Parameters:
filename_out : str

Path to where the results should be saved.

shuffle : bool, optional

If True, events are shuffled before being saved. That’s important when there are multiple distinct samples (e.g. signal and background). Default value: True.

Returns:
None
set_acceptance(self, pt_min_e=None, pt_min_mu=None, pt_min_a=None, pt_min_j=None, eta_max_e=None, eta_max_mu=None, eta_max_a=None, eta_max_j=None)[source]

Sets acceptance cuts for all visible particles. These are taken into account before observables and cuts are calculated.

Parameters:
pt_min_e : float or None, optional

Minimum electron transverse momentum in GeV. None means no acceptance cut. Default value: None.

pt_min_mu : float or None, optional

Minimum muon transverse momentum in GeV. None means no acceptance cut. Default value: None.

pt_min_a : float or None, optional

Minimum photon transverse momentum in GeV. None means no acceptance cut. Default value: None.

pt_min_j : float or None, optional

Minimum jet transverse momentum in GeV. None means no acceptance cut. Default value: None.

eta_max_e : float or None, optional

Maximum absolute electron pseudorapidity. None means no acceptance cut. Default value: None.

eta_max_mu : float or None, optional

Maximum absolute muon pseudorapidity. None means no acceptance cut. Default value: None.

eta_max_a : float or None, optional

Maximum absolute photon pseudorapidity. None means no acceptance cut. Default value: None.

eta_max_j : float or None, optional

Maximum absolute jet pseudorapidity. None means no acceptance cut. Default value: None.

Returns:
None