SpartyJet documentation
hosted by CEDAR HepForge
Last updated: Wed Mar 10 20:06:15 2010

Overview

Spartyjet is a set of software tools to run jet finding algorithms.
Besides physics motivations, it has been written with 3 goals in mind :

  • ease of use
  • flexibility
  • speed
To meet these goals, we tried to write it in a modular way; it consists of bricks of software put together allowing one:
  • to include any kind of input in the form of 4-vectors
  • to perform any operation on this 4-vector set (including jet finding of course)
  • to save all results and any associated quantities into a ROOT ntuple
See a schematic view of the code organization here.
A lots of different algorithm are shipped with Spartyjet : Atlas, CMS, CDF and D0 jet algorithms, the FastJet library, and Pythia's CellJet.

Installation

Requirements

OS
Most linux distributions, Mac OSX.
Compiler
GCC/g++ (with make).
ROOT
Recent version, 5.18 or later (see notes).
python
python2.4 or later (not entirely necessary but highly recommended)
fortran compiler (optional)
Needed to compile the libraries allowing StdHep reading.

Compilation

  • In spartyjet/ directory, just type:
    make
    

  • Libraries
    This will build a list of libraries in spartyjet/libs that you can load from a ROOT session or Python script, or you can link to from an executable.

    libs/libJetCore.so - Core infrastructure
    libs/libJetTools.so - JetTools
    libs/libCDFJet.so - CDF JetClustering and MidPoint
    libs/libD0Jet.so - D0 cone algortihms
    libs/libATLASJet.so - ATLAS algorithms
    libs/libCellJet.so - CellJet algorithm
    libs/libFastJet.so - FastJet and SISCone
    libs/libEventShape.so - Thrust and other eventshapes
    libs/libSpartyDisplay.so - SpartyJet GUI

  • Options: SpartyJet has two buidling options to note:
    1. StdHEP libraries: These are optional to compile. If you plan on using StdHEP files as input (and you have a fortran compiler), then open spartyjet/MakefileCommonOpt and uncomment the line "DO_STDHEP =1" and specify which fortran you want to use with "F77 = ". On the other hand, if StdHEP compilation is producing errors, then make sure DO_STDHEP is commented out. The current default is to compile StdHEP.
    2. Fastjet FastJet: SpartyJet ships with the latest version of FastJet. If however, you prefer to use your own installation, you can pass SpartyJet the location of the FastJet installation. SpartyJet will then link to this rather than compiling the default version. Simply open spartyjet/MakefileCommonOpt and specify the install location with "FASTJETDIR = /path/to/fastjet/installation". NOTE: If you have linking problems between your version of fastjet and SpartyJet, either recompile your fastjet with the "-with-pic" option enabled before compilation, or allow SpartyJet to compile its own version.

Running

  • Working examples of how to use SpartyJet, can be found in the following directories:
    spartyjet/examples_python : Python scripts RECOMMENDED
    spartyjet/examples_ROOT : ROOT scripts
    spartyjet/examples_C : Compiled C++ programs

  • The python interface to SpartyJet is the most powerful. To use it, please execute:
    source setup.sh
    
    in the spartyjet/ directory. This allows the spartyjet python modules and libraries to be accessible from any directory.

  • All examples expect to be run from their own directory as follows:
    Python script
    cd examples_py
    python simpleExample.py
    
    ROOT script
    cd examples_ROOT
    root -l simpleExample.C
    
    Compiled C++
    cd examples_C
    ./simpleExample.exe
    

  • Detailed walk-throughs of simple examples can be found in the Examples Explained section

Input

All SpartyJet jobs need an InputMaker object to read input from some data file and prepare a list of 4-vectors for JetTools to process. There are several types of InputMakers available. An example of the implementation for each type of input can be seen in:

  • examples_py/inputExample.py
  • examples_ROOT/inputExample.cc
  • examples_C/inputExample.C

NtupleInputMaker

Sample: data/data/J2_clusters.root

This form of input reads ROOT files. This InputMaker requires the components of the input 4-vectors to be stored in separate branches of a ROOT TTree. The following definitions are supported:

  • px,py,pz,e
  • (psuedo)rapidity,phi,pt,e
  • (psuedo)rapidity,phi,pt,m
You need only specify how the information is stored (array or vector / float or double) and the names of the branches.

As an example, to set up an NtupleInputMaker to read the file spartyjet/data/J2_clusters.root If you don't know how variables are internally stored in your ntuple, do the following:

 
root -l data/J2_clusters.root
root [0] clusterTree->MakeClass("test");
Open the file test.h and check to see how the variables are stored. In this example, we see lines like:
 
vector<float>   *Cluster_eta;
indicating that our 4-vectors are stored in vectors of floats. Now to configure SpartyJet to accept this, we need to:

  • Create an NtupleInputMaker of the correct type: (in our case: vector, float, (eta,phi,pt,E) )
    input = SJ.NtupleInputMaker(SJ.NtupleInputMaker.EtaPhiPtE_vector_float)
    
    For full list of Input codes see: JetCore/InputMaker_Ntuple.hh

  • Configure the names of the TBranches
    input.set_prefix('Cluster_')
    input.set_n_name('N')
    input.set_variables('eta','phi','p_T','e')
    input.setFileTree('../data/J2_clusters.root', 'clusterTree')
    

  • Specify if input is massless, only useable in eta,phi,ptm,E mode (if true, pt is ignored):
    input.set_masslessMode(True)
    

  • Set the input file and tree names:
    input.setFileTree('../data/J2_clusters.root', 'clusterTree')
    

Python Shortcut:
To allow SpartyJet to configure your ntuple using some assumptions use the following helper function:

input = createNtupleInputMaker('../data/J2\_clusters.root', inputprefix='Cluster')

StdTextInput

Sample: data/J1_Clusters.dat

This form of input reads ASCII files.
To separate events, put one of the following lines between the events

.Event .event (only the .E or .e is important) N n

The form of the four vectors should be:

E px py pz

This input is configured simply with:

input = SJ.StdTextInput('../data/J1_Clusters.dat')

If the form is the opposite (px py pz E), then call the function

input.invert_input_order(True);
and it will be read in properly.

An example of this input can be seen in data/J1_Clusters.dat

StdHepInput

Sample: data/ttbar_smallrun_pythia_events.hep

This form of input reads StdHEP format XDR files. It will look for particles with the status code of 1 (final state). It is also able extract the pdgId from the input data to allow further filtering and matching.

This input is configured simply with:

SJ.StdHepInput('../data/ttbar_smallrun_pythia_events.hep')

NOTE: To read StdHep, one must turn on the StdHep compilation as explained here.

CalchepPartonTextInput

Sample: data/gg_ggg_events.dat

This form of input reads output from CALCHEP. It reads in the number of initial and final state particles, and then for each event saves only the information for the final state particles.

This input is configured simply with:

SJ.CalchepPartonTextInput('../data/gg_ggg_events.dat')

HepMCInput

Sample: data/HepMC_sample.dat

This form of input reads HepMC format ascii files. HepMC ascii files contain the following separators:

E - Denotes new event
V - Information about a vertex
P - Information about a particle

This class reads in the four vectors of the particles denoted with a status code of 1 (not decayed, final state).


Input Options

Rejecting bad input
The InputMaker can be set to remove 4-vectors with negative energy and non-physical momenta by using the following function.

input.reject_bad_input(bool)
The current default is false; no checks will be done. (An alternative to this method of dealing with bad input is given in the JetBuilder Options section).

Reading of ParticleDataGroup Id Codes
The InputMaker can be set to read pdg Id codes from the input data with the following function.

input.readPdgId(bool)
This makes the pdg Id available for input selection as well as persistifying the pdg Id's of the input particles for offline ananlysis.


JetBuilder

JetBuilder is the job manager for SpartyJet. As such, it takes the input from one of the above InputMakers and passes it through a sequence of JetTools, including any JetAlgorithms. It then gives the final list of jets to NtupleMaker which prepares the output.

Input Functions

Input is passed to JetBuilder via:
configure_input(input)
where input is an InputMaker object.

Minimum Bias Overlay

Example: examples_py/overlayExample.py

JetBuilder allows the user to add minimum bias events to the signal events to study the effects of pileup. To enable minimum bias events, one must:

  1. Create another input object of any type (StdTextInput,NtupleInput...):
    MBinput = SJ.StdTextInput('data/MB_Clusters.dat')
    
  2. Tell JetBuilder to add n minbias events to each signal event:
    builder.add_minbias_events(n,MBinput);
    
JetBuilder will start at the beginning of the minbias data file and read the first n events for the first data event, then the second n events for the second data event, so on. When the end of the minbias file is found, it will simply continue from the beginning.

Note: There is a third optional boolean argument to add_minbias_events() that tells JetBuilder whether to draw the number of minbias events from a Poission distribution. In this case, n is the Possion mean or expected number of MinBias events.

JetTool Options

JetTools are blocks of code that act on a set of Jets. They are defined by the user and then passed either to JetBuilder (this executes the JetTool on all algorithms) or to a specific Jet Algorithm, where they will be executed only on that particular algorithm's jets. Details on specific tools can be found in the JetTools section.

JetTools are added to the exection sequence via:

add_jetTool()
or to the front of the sequence with:
add_jetTool_front()

JetAlgorithms are added to the via:

add_default_alg(alg)
This adds the jet algorithm alg to the job and runs on its jets all the tools that have been passed to the JetBuilder.

Alternatively, to add a JetAlgorithm with only the tools explicitly associated to it use:

add_custom_alg(alg)

Constituents

SpartyJet has the ability to retain the information about each jet's constituents. The method by which this information is stored is explained in the output section By default, this is enabled, to disable constituent saving, you must add false as a second argument when using builder.add_default_alg()

Note: The CellJet algorithm is currently incapable of determining constituent information, so if constituents are turned on, the branches in the ROOT TTree output will be empty.

Text File Output

The text file output option tells JetBuilder to produce a text file that contains a list of all jets found from all algorithms for all events.
The file is easy to read, and recommended for quick viewing of results.
To turn on the text file output, you must call builder.add_text_output() , and pass it the filename you want to create. For example:

builder.add_text_output("data/output/text_output.dat")

Rejecting Bad Input

To deal with bad input 4-vectors, JetBuilder makes use of a JetTool called NegEnergyTool. It can call this function before and after the jet finder tool to:

  1. Store which constituents have negative energy
  2. Inverse their sign for jet finding
  3. Check if any constiuents of each jet had negative energy and correct accordingly
The default is to NOT use this tool, and one may switch it on with: builder.do_correct_neg_energy(bool)


JetTools

Jet Tools are blocks of code that act on a JetCollection every event. They generally before one or more of the following function:
  • Add or remove jets from the list
  • Modify the jets themselves
  • Add information about each jet as a Jet Moment
  • Add information about each event as an Event Moment
Jet tools can be added either before or after the primary Jet finder (ex Kt algorithm) in the job. The following is a list of Jet Tools shipped with SpartyJet.

Selectors

These tools allow the user to remove some Jets from an input/output JetCollection:

  • JetPtSelectorTool(double ptmin) - Removes jets with $p_{T}$ below cut
  • JetPtORESelectorTool(doule ptmin,double emin) - Removes jets with $p_{T}$ or energy below cuts
  • JetEtaCentralSelectorTool(double abs_eta) - Removes jets outside of central eta region
  • JetEtaForwardSelectorTool(double abs_eta) - Removes jets inside of central eta region
  • JetInputPdgIdSelectorTool(std::vector<int> pdgIds) - Removes input jets with pdgIds

HullMomentTool

This tool finds the convex hull enclosing each jet and saves the hull length and area as "Hull" and "HullA" respectively.

EtaPhiMomentTool

This tool calculates angular 2nd moments in eta and phi of each jet stores them as "M2eta" and "M2phi"

PtDensityTool

This tool calculates each event's $p_{T}$ density using fastjet. It does this by finding all the jets in the event with no min $p_{T}$ requirement. And then extracting the ptDensity from thse jets by selection the median pTDensity for each bin in $eta$. These ptDensities are stored as "ptDensity" and the eta bin limits are stored as "ptDensityBins"

JetAreaCorrectionTool

This tool uses the area of each jet and the pTDensity found by the PtDensityTool to calcualte a correction (stored as jet moment JetAreaCorr) to the Jet's $p_{T}$.


YSplitterTool

This tool uses fastjet to calculate the yvalues asociated with a set of recombinations. In this implementation, SpartyJet will run a fastjet algorithm of the User's choice on the constituents of a given jet. It can be called with either of the following constructors:
YSplitterTool(float R, fastjet::JetAlgorithm alg, int ny, int njet) YSplitterTool(fastjet::JetDefinition *jet_def, int ny, int njet) where n_jets is the number of jets for which the y_values will be calculated and n_y_values is the number of y_values to calculate for each jet.

EConversionTool

This tool simply converts the units of all the jets between MeV and GeV. The user can convert to arbitrary units as well.
EConversionTool("ConversionTool",bool toGeV)

Algorithms

Here is a list of the algorithms with a list of the options and parameters that can be set.
All of these Algorithms can be added as explained in JetBuilder.

Fastjet-specific features

Output

When using the JetBuilder class the result is a ROOT ntuple containing jet variables for each algorithm added, plus variables for input particles to the jet algorithms.

Output variable type

It is possible to choose the type of the variable saved in the ntuple build by JetBuilder. Choice are between pure C array or STL vector of float or double.

builder.output_var_style.array_type ="array" // (default) other option is "vector"
builder.output_var_style.array_type ="float" // (default) other option is "double"

Constituents

For all algorithms, (except CellJet in v1), constituents information can be saved as follows.
Assuming a jet finder named "MyJet" has been added, one gets 2 variables in the ROOT TTree:
MyJet_numC
MyJet_ind

MyJet_numC is an array of size MyJet_N
MyJet_numC[i] is the number of constituents of ith jet.

MyJet_ind is an array of size InputJet_N .
MyJet_ind[i] is the index of the jet to which the ith input constituent has been assigned.

For example in ROOT:

myTree->Draw("InputJet_e","FastKt_ind==0")
will give the energy distribution of constituents in jet number 0 (i.e. highest pt jets) for the FastKt collection.

Analysis

There is currently work being done on a graphical user interface for intereactive SpartyJet analysis.
Users can test this out with the guiExample.py in the examples_py directory.


Examples Explained

The following walkthrough follows the script: spartyjet/python_py/simpleExample.py

Goal: Run the FastJet AntiKt Algorithm on the first 10 events listed in data/J1_Clusters.dat

  • Load the libraries that are needed for the algorithms that you are running.
    from ROOT import *
    gSystem.Load("libPhysics.so")
    gSystem.Load("../libs/libJetCore.so") # All jobs need this
    gSystem.Load("../libs/libFastJet.so") # Needed for AntiKt alg
    from ROOT import SpartyJet as SJ
    

    In practice for python scripts, all necessary includes can be taken care of with the following line: (You need to source setup.sh in the spartyjet directory first)

    from  SpartyJetConfig import *
    

  • Create JetBuilder object to manage SpartyJet job. The argument sets the output message level (DEBUG,INFO,WARNING,ERROR).
    builder = SJ.JetBuilder(SJ.INFO)
    

  • Create an input object of type StdTextInput with the filename containing the events.
    SJ.StdTextInput("data/J1_Clusters.dat");
    

  • Pass input object to JetBuilder.
    builder.configure_input(input)
    

  • Define Jet algorithm that you want to run.
    antikt4 = SJ.FastJet.FastJetFinder('AntiKt4',fj.antikt_algorithm,0.4)
    

  • Pass algorithm to your JetBuilder
    builder.add_default_alg( antikt4 )
    

  • Configure (optional) simple text output for quick visual check of results.
    builder.add_text_output("data/output/text_simple.dat");
    

  • Configure the ntuple output by specifying the name of the tree and the ROOT file you want the data to be stored in.
    builder.configure_output("SpartyJet_Tree","data/output/simple.root");
    

  • Give the command to run the algorithms on the first 10 events
    builder.process_events(10);
    

  • Then running the script will process the first 10 events on the file specified in the input object and produce the .root file specified in configure_output. This is meant as a basic introduction, and there are a lot more functions than listed here. Look at the other examples particularly in the python directory.

Notes

MidPoint second pass:
If you have an error related to midpoint second pass, then call

MidPointFinderObject.set_SavePassnum(false);
This will turn off the saving of the pass number information, which has been seen to cause problems in some cases

Bad Input:
There are two ways to deal with bad input to SpartyJet, one is described in "Section" and one is described in the "JetBuilder section"

Any questions/comments/suggestions, email:

Pierre-Antoine Delsart: delsart@dontspamin2p3.fr
Joey Huston: huston@dontspam.pa.msu.edu
Brian Martin: marti347@msudontspam.edu

About this document ...

SpartyJet documentation

This document was generated using the LaTeX2HTML translator Version 2008 (1.71)

Copyright © 1993, 1994, 1995, 1996, Nikos Drakos, Computer Based Learning Unit, University of Leeds.
Copyright © 1997, 1998, 1999, Ross Moore, Mathematics Department, Macquarie University, Sydney.

The command line arguments were:
latex2html -split 0 -dir html -no_navigation SpartyJet_Documentation.tex

The translation was initiated by Brian Martin on 2010-03-10


Brian Martin 2010-03-10