Skip to main content

Documentation Index

Fetch the complete documentation index at: https://mintlify.com/ShipSoft/FairShip/llms.txt

Use this file to discover all available pages before exploring further.

macro/ShipAna.py is the standard post-reconstruction analysis script for FairShip. It opens both the MC simulation file and the matching reconstruction file, joins them via ROOT friend trees, and loops over events to fill histograms characterising track quality, vertex resolution, momentum pulls, and HNL candidate kinematics. A complementary, higher-level interface is provided by the experimental analysis_toolkit module demonstrated in examples/analysis_example.py. This page documents both paths — the low-level ShipAna.py approach and the toolkit approach — so you can choose the one that best fits your workflow.

Running ShipAna

Command-line reference

-f / --inputFile
string
required
Path to the MC simulation ROOT file containing the cbmsim tree. Accepts a comma-separated list of files, which are opened as a TChain.
-g / --geoFile
string
required
Path to the geometry ROOT file. Used to initialise detector modules, the magnetic field map, and GenFit material effects — all needed for track extrapolation.
-r / --recoFile
string
Path to the reconstruction file (*_rec.root). If omitted, ShipAna.py automatically derives the name by replacing .root with _rec.root in the input filename. Accepts a comma-separated list of files when -f is also a list.
-n / --nEvents
integer
default:"999999"
Maximum number of events to analyse. The actual count is min(sTree.GetEntries(), nEvents).
--Debug
flag
Enable verbose per-event debug printing.

Usage examples

python macro/ShipAna.py \
    -f sim_my-run.root \
    -g geo_my-run.root

Friend-tree setup inside ShipAna

ShipAna.py establishes the MC ↔ reco connection in the file-opening block. For a single file:
f = ROOT.TFile(options.inputFile)
sTree = f["cbmsim"]
sTree.AddFriend("ship_reco_sim", options.recoFile)
For comma-separated multi-file inputs:
sTree    = ROOT.TChain("cbmsim")
recoChain = ROOT.TChain("ship_reco_sim")
for x in options.inputFile.split(","):
    sTree.AddFile(x)
for x in options.recoFile.split(","):
    recoChain.AddFile(x)
sTree.AddFriend(recoChain)
After AddFriend, every branch in ship_reco_sim is accessible directly on sTree. The event loop then calls sTree.GetEntry(n) once per event and all branches — MC and reco — are loaded together.

The cbmsim tree structure

ShipAna.py reads the following branches from the combined cbmsim + ship_reco_sim view:
BranchSource treeDescription
MCTrackcbmsimMonte Carlo particle list; index 0 is the primary, index 1+ are daughters
strawtubesPointcbmsimTrue MC hit points in the straw spectrometer
ShipEventHeadership_reco_simRun/event metadata including event time
FitTracksship_reco_simgenfit::Track* fitted track objects
fitTrack2MCship_reco_simInteger vector: reco track index → MC track index
Particlesship_reco_simShipParticle decay-vertex candidates
Trackletsship_reco_simShort track segments used during fitting
When the pattern recognition was run with the _PR variant (from an earlier workflow), ShipAna automatically redirects:
if sTree.GetBranch("FitTracks_PR"):
    sTree.FitTracks = sTree.FitTracks_PR
if sTree.GetBranch("fitTrack2MC_PR"):
    sTree.fitTrack2MC = sTree.fitTrack2MC_PR
if sTree.GetBranch("Particles_PR"):
    sTree.Particles = sTree.Particles_PR

Quality cuts

Four module-level constants in ShipAna.py control which tracks and vertices are accepted:
VariableDefaultMeaning
chi2CutOff4.0Maximum χ²/NDF for a track to enter momentum-resolution histograms
fiducialCutFalseIf True, require decay vertex inside the fiducial decay volume
measCutFK25Minimum number of GenFit measurement degrees of freedom (full Kalman)
measCutPR22Minimum NDF when using the _PR pattern-recognition tracks
docaCut2.0Maximum DOCA (cm) between the two tracks of a vertex candidate
measCut is set to measCutFK by default and switched to measCutPR automatically if the FitTracks_PR branch is found. The track loop only enters the chi2 and IP histograms for tracks that pass the NDF cut; only vertices with doca < docaCut enter the HNL mass histogram.

Histograms produced

ShipAna.py books and fills the following histogram families:
KeyDescription
delPOverP(p_truth − p_reco) / p_truth vs p_truth
delPOverPz(1/pz_truth − 1/pz_reco) × pz_truth vs pz_truth
delPOverP2Same as delPOverP but only for χ²/NDF < chi2CutOff
pullPOverPx/y/zMomentum pulls per component vs p_truth

Output files

At the end of the event loop, makePlots() renders canvas plots and ut.writeHists() writes all histograms. The output filename is derived from the input:
sim_my-run.root       →  sim_my-run_ana.root
sim_my-run_rec.root   →  sim_my-run_ana.root   (strips _rec before adding _ana)
When the input path begins with /eos/ or is a comma-separated list, the file is written to the current working directory instead of alongside the input.

The experimental analysis_toolkit

python/experimental/analysis_toolkit.py provides a higher-level object-oriented interface for signal selection studies. It is demonstrated in examples/analysis_example.py.

Opening files and adding the friend tree

import ROOT
from experimental import analysis_toolkit

f = ROOT.TFile.Open("sim_my-run.root", "read")
tree = f["cbmsim"]
tree.AddFriend("ship_reco_sim", "sim_my-run_rec.root")

geo_file = ROOT.TFile.Open("geo_my-run.root", "read")

selection_check

analysis_toolkit.selection_check(geo_file) initialises geometry-aware selection helpers from the geo ROOT file. It loads FAIRGeom and the ShipGeo dictionary, reads the veto geometry YAML for the configured decay-volume medium (helium or vacuum), and provides the following methods on each ShipParticle candidate:
MethodReturnsDescription
access_event(tree)NoneMust be called once per event before any per-candidate methods
define_candidate_time(candidate)float (ns)Reconstructs decay-vertex time from straw MC hits and propagation speed
impact_parameter(candidate)float (cm)IP of candidate 4-momentum relative to (0, 0, target.z0)
dist_to_innerwall(candidate)float (cm)Distance from vertex to inner wall of decay vessel
dist_to_vesselentrance(candidate)float (cm)Distance from vertex to upstream entrance of decay vessel
DOCA(candidate)float (cm)Distance of closest approach between the two daughter tracks
invariant_mass(candidate)float (GeV)Di-track invariant mass
daughtermomentum(candidate)(float, float) (GeV)Momenta of daughter 1 and daughter 2
nDOF(candidate)(int, int)NDF of daughter track fits
chi2nDOF(candidate)(float, float)χ²/NDF of daughter track fits
preselection_cut(candidate, IP_cut, show_table)boolCombined pre-selection flag

event_inspector

analysis_toolkit.event_inspector() provides a lightweight event-dumping utility:
inspector = analysis_toolkit.event_inspector()
inspector.dump_event(event, mom_threshold=0.5)  # prints tracks with p > 0.5 GeV

Complete analysis_example workflow

#!/usr/bin/env python
# examples/analysis_example.py
import ROOT
import rootUtils as ut
from experimental import analysis_toolkit

f = ROOT.TFile.Open("sim_my-run.root", "read")
tree = f["cbmsim"]
tree.AddFriend("ship_reco_sim", "sim_my-run_rec.root")

geo_file  = ROOT.TFile.Open("geo_my-run.root", "read")
selection = analysis_toolkit.selection_check(geo_file)
inspector = analysis_toolkit.event_inspector()

hist_dict = {}
ut.bookHist(hist_dict, "inv_mass",        "Invariant mass; GeV",      200, 100, 100)
ut.bookHist(hist_dict, "impact_parameter","Impact parameter; cm",      200, 100, 100)
ut.bookHist(hist_dict, "DOCA",            "DOCA; cm",                  200, 100, 100)

for event_nr, event in enumerate(tree):
    selection.access_event(event)
    if len(event.Particles) == 0:
        continue

    inspector.dump_event(event, mom_threshold=0.5)

    for candidate_id, signal in enumerate(event.Particles):
        hist_dict["inv_mass"].Fill(selection.invariant_mass(signal))
        hist_dict["impact_parameter"].Fill(selection.impact_parameter(signal))
        hist_dict["DOCA"].Fill(selection.DOCA(signal))

        IP_cut = 250
        if selection.preselection_cut(signal, IP_cut=IP_cut, show_table=False):
            print(f"Event {event_nr}, candidate {candidate_id}: passes pre-selection")

ut.writeHists(hist_dict, "preselectionparameters.root")
The analysis_toolkit reads the decay-volume medium from ShipGeo.DecayVolumeMedium and automatically loads the correct veto geometry YAML (veto_config_helium.yaml or veto_config_vacuums.yaml). Make sure the FAIRSHIP environment variable points to the FairShip source directory so the YAML can be found.
The analysis_toolkit module lives under python/experimental/ and its API may change between releases. Pin the FairShip commit hash in any analysis that relies on specific method signatures.

Build docs developers (and LLMs) love