FunTuple.functorcollections

Functor collections for the FunTuple framework.

Refer to the docstrings of the collections in this submodule for examples of typical usage.

Functions

EventInfo(→ FunTuple.FunctorCollection.FunctorCollection)

Event-level collection of most-used functors for event information.

SelectionInfo(...)

Event-level collection of most-used functors for tupling trigger/Sprucing information.

HltTisTos(→ FunTuple.FunctorCollection.FunctorCollection)

Candidate-level collection to store TIS (Trigger Independent of Signal) or TOS (Trigger On Signal)

Kinematics(→ FunTuple.FunctorCollection.FunctorCollection)

Candidate-level collection of most-used functors on kinematics.

MCKinematics(...)

Candidate-level collection of most-used functors on kinematics on the related MC particle.

DecayTreeFitterResults(...)

Candidate-level collection of functors that access DecayTreeFitter fit results.

MCPromptDecay(...)

Candidate-level collection of functors to determine if a decay is prompt based on the true lifetime of its ancestors.

MCHierarchy(→ FunTuple.FunctorCollection.FunctorCollection)

Candidate-level collection of functors yielding the hierarchy/relationships among MC particles.

MCVertexInfo(...)

Candidate-level collection of functors on vertices of the related MC particle.

MCPrimaries(→ FunTuple.FunctorCollection.FunctorCollection)

Event-level collection of functors to store the properties (X,Y,Z,T) of all MC primary vertices.

MCReconstructible(...)

Candidate-level collection of functors to store reconstrutible information for the MC particle.

MCReconstructed(...)

Candidate-level collection of functors to store reconstrutible information for the MC particle.

ParticleIsolation(...)

Candidate-level collection of functors on track or neutral isolation, using information from relations between sets of particles determined by a given criterion.

ConeIsolation(...)

Candidate-level collection of functors on charged and neutral isolation, using information from relations between sets of particles determined by a given criterion.

VertexIsolation(...)

Candidate-level collection of functors on vertex isolation, using information from relations between sets of particles determined by a given criterion.

NeutralCaloInfo(...)

Candidate-level collection of functors on neutral calorimeter hypotheses information.

ChargedCaloInfo(...)

Candidate-level collection of functors on charged calorimeter hypotheses information.

ParticleID(→ FunTuple.FunctorCollection.FunctorCollection)

Candidate-level collection of functors on global particle identification variables.

RecSummary(→ FunTuple.FunctorCollection.FunctorCollection)

Event-level collection providing access to the RecSummary record.

FlavourTaggingResults(...)

Candidate-level collection of functors that access FlavourTagging results.

Module Contents

EventInfo() FunTuple.FunctorCollection.FunctorCollection[source]

Event-level collection of most-used functors for event information. By default the FunTuple framework stores the RUN_NUMBER and EVENT_NUMBER variables. When running over truth samples (MCParticles) i.e. xgen files, the RUN_NUMBER and EVENT_NUMBER are not yet available (see issue https://gitlab.cern.ch/lhcb/DaVinci/-/issues/90).

Output variables:
  • BUNCHCROSSING_ID : F.BUNCHCROSSING_ID(odin)

  • BUNCHCROSSING_TYPE : F.BUNCHCROSSING_TYPE(odin)

  • ODINTCK : F.ODINTCK(odin)

  • GPSTIME : F.GPSTIME(odin)

Example:

from FunTuple import FunTuple_Particles as Funtuple
import FunTuple.functorcollections as FC
from DaVinci import options

variables = FC.EventInfo()
# ...

tuple = Funtuple(name="MyTuple",
                 fields=...,
                 variables=variables,
                 ...
                 )
SelectionInfo(*, selection_type: GaudiConf.LbExec.HltSourceID, trigger_lines: list[str]) FunTuple.FunctorCollection.FunctorCollection[source]

Event-level collection of most-used functors for tupling trigger/Sprucing information.

Parameters:
  • selection_type (HltSourceID) – Name of the selection type i.e. “Hlt1” or “Hlt2” or “Spruce”. Used as branch name prefix when tupling and as source ID to get decision reports.

  • trigger_lines (list(str)) – List of line names for which the decision is requested.

Output variables:
  • _TCK : F.TCK(dec_report) - for each source ID

Example:

from FunTuple import FunTuple_Particles as Funtuple
from DaVinci import options
import FunTuple.functorcollections as FC

# List of lines you want decisions for
selection_type = HltSourceiD.Hlt2
hlt2_lines = ['Hlt2CharmD0ToKmPipLineDecision']

variables = FC.SelectionInfo(selection_type, hlt2_lines)

tuple = Funtuple(name="MyTuple",
                 fields=...,
                 variables=variables,
                 ...
                 )
HltTisTos(*, selection_type: GaudiConf.LbExec.HltSourceID, trigger_lines: list[str], data: PyConf.dataflow.DataHandle, trigger_particles: list[PyConf.dataflow.DataHandle] | None = None) FunTuple.FunctorCollection.FunctorCollection[source]

Candidate-level collection to store TIS (Trigger Independent of Signal) or TOS (Trigger On Signal) information of candidates for specified Hlt1/Hlt2 lines. For TOB (Trigger On Both), the condition is (!TIS && !TOS && event_decision). To store “event decision” of a line use the SelectionInfo functor collection.

Parameters:
  • selection_type (HltSourceID) – Name of the selection type i.e. “Hlt1” or “Hlt2”. Currently only “Hlt1” is supported. Used as branch name prefix when tupling and as source ID to get decision reports.

  • trigger_lines (list(str)) – List of Hlt1/Hlt2 line names for which the decision is requested.

  • data (DataHandle) – TES location of particles from a Sprucing/Hlt2 line i.e. output of for example get_particles("/Event/HLT2/<linename>/Particles").

  • trigger_particles (list(DataHandle), optional) – List of TES location of trigger particles. This could be specified for example when running over “Hlt2” process directly (where only TOS information might be useful). When running over “Spruce” or “TurboPass” process this is not needed as the path is inferred from the DV configuration.

Output variables:
  • _TOS : F.IS_TOS(line, TisTosRelations) - for each provided line

  • _TIS : F.IS_TIS(line, TisTosRelations) - for each provided line

Example:

from DaVinci import options
import FunTuple.functorcollections as FC

# Get input data from a particular "<linename>" (can be Sprucing line or Hlt2 line)
data = get_particles(f"/Event/HLT2/<linename>/Particles")

# List of lines for which Tis/Tos information is requested
trigger_lines = ["Hlt1TwoTrackMVADecision", "Hlt1TrackMVADecision"]
variables = FC.HltTisTos(selection_type="Hlt1", trigger_lines=trigger_lines, data=data)

# Store event decision to check if a candidate is "TOB" wrt to a given line
evt_variables = FC.SelectionInfo(selection_type="Hlt1", trigger_lines)

tuple = Funtuple(name="MyTuple",
                 fields=...,
                 variables=variables,
                 event_variables=evt_variables,
                 ...
                 )
Kinematics() FunTuple.FunctorCollection.FunctorCollection[source]

Candidate-level collection of most-used functors on kinematics.

Output variables:
  • M : F.MASS

  • P : F.P

  • P[T,X,Y,Z] : F.P[T,X,Y,Z]

  • ENERGY : F.ENERGY

Example:

from FunTuple import FunTuple_Particles as Funtuple
import FunTuple.functorcollections as FC

variables = FC.Kinematics()

tuple = Funtuple(name="MyTuple",
                 fields=...,
                 variables=variables,
                 ...
                 )
MCKinematics(*, mctruth_alg: DaVinciMCTools.MCTruthAndBkgCat) FunTuple.FunctorCollection.FunctorCollection[source]

Candidate-level collection of most-used functors on kinematics on the related MC particle.

Parameters:

mctruth_alg (DaVinciMCTools.MCTruthAndBkgCat) – An instance of the MCTruthAndBkgCat helper class that is used for retrieving MC truth information of a reconstructed particle.

Output variables:
  • TRUEP : mctruth_alg(F.P)

  • TRUEPT : mctruth_alg(F.PT)

  • TRUEP[X,Y,Z] : mctruth_alg(F.P[X,Y,Z])

  • TRUEENERGY : mctruth_alg(F.ENERGY)

Example:

from FunTuple import FunTuple_Particles as Funtuple
import FunTuple.functorcollections as FC
from DaVinciMCTools import MCTruthAndBkgCat

# Load particles onto TES location for a given event cycle
input_data = get_particles("/Event/HLT2/<linename>/Particles")

# Define an algorithm that builds one-to-one relation table b/w Reco Particle -> Truth MC Particle.
# Here "options" is the DV options object (see DaVinciExamples or DaVinciTutorials)
mctruth_alg = MCTruthAndBkgCat(input_data)

variables = FC.MCKinematics(mctruth_alg)

tuple = Funtuple(name="MyTuple",
                 fields=...,
                 variables=variables,
...
)
DecayTreeFitterResults(*, DTF: DecayTreeFitter.DecayTreeFitter, decay_origin: bool = True, with_lifetime: bool = False, with_kinematics: bool = True, prefix: str = 'DTF') FunTuple.FunctorCollection.FunctorCollection[source]

Candidate-level collection of functors that access DecayTreeFitter fit results.

Parameters:
  • DTF (DecayTreeFitter.DecayTreeFitter) – An instance of the DecayTreeFitter helper class that builds the necessary functors.

  • decay_origin (bool) – Returns functorcollection for decay origin. Defaults to True.

  • with_lifetime (bool, optional) – The returned functors comprise lifetime and flight distance. Defaults to False.

  • with_kinematics (bool, optional) – The returned functors comprise kinematics. Defaults to True.

  • prefix (str, optional) – Additional prefix to output functors. Defaults to ‘DTF’.

Output variables:
  • _CHI2

  • _NDOF

  • _CHI2DOF

  • _NITER

  • _MASS

  • _MASSERR

  • _P

  • _PERR

If with PV constraint:
  • _PV_KEY

  • _PV_[X,Y,Z]

If with lifetime:
  • _CTAU

  • _CTAUERR

  • _FD

  • _FDERR

If with kinematics:
  • _P[X,Y,Z,E]

Example:

from FunTuple import FunTuple_Particles as Funtuple
import FunTuple.functorcollections as FC
from DecayTreeFitter import DecayTreeFitter

DTF = DecayTreeFitter(
    name = 'DTF',
    input_particles = get_particles(...),
    input_pvs = get_pvs()
)

variable_decay_origin = FC.DecayTreeFitterResults(DTF, decay_origin=True, with_lifetime=True, with_kinematics=True)
variable_kinematics = FC.DecayTreeFitterResults(DTF, decay_origin=False, with_lifetime=False, with_kinematics=True)

variables = {
    'B': variable_decay_origin,
    'phi_1020': variable_kinematics,
    'Kplus': variable_kinematics,
    'Kminus': variable_kinematics,
}

tuple = Funtuple(name="MyTuple",
                 fields=...,
                 variables=variables,
                 ...
                 )
MCPromptDecay(*, mctruth_alg: DaVinciMCTools.MCTruthAndBkgCat | None = None, max_lifetime: float = 1e-07 * SystemOfUnits.ns, store_longlived_particle_id: bool = True) FunTuple.FunctorCollection.FunctorCollection[source]

Candidate-level collection of functors to determine if a decay is prompt based on the true lifetime of its ancestors. If a decay does not have any long-lived ancestors, it is classified as prompt.

This algorithm is based on the Runs 1 and 2 MCTupleToolPrompt tool (https://gitlab.cern.ch/lhcb/Analysis/-/blob/run2-patches/Phys/DecayTreeTupleMC/src/MCTupleToolPrompt.cpp).

If no argument is provided, i.e. mctruth_alg==None, this collection creates functors that apply to MCParticles; otherwise, the functors created are assumed to apply to Particles, which will find the associated MCParticle automatically using the input truth matching.

Parameters:
  • mctruth_alg (DaVinciMCTools.MCTruthAndBkgCat, optional) – An instance of the MCTruthAndBkgCat helper class that is used for retrieving MC truth information of a reconstructed particle. Defaults to None.

  • max_lifetime (float, optional) – maximum lifetime of short-lived particles (ns). Defaults to 1e-7 ns.

  • store_longlived_particle_id (bool, optional) – create functors to store the particle ID and container key of the first long-lived ancestor. Defaults to True.

Output variables:
  • MC_ISPROMPT

  • MC_LONGLIVED_ID

  • MC_LONGLIVED_KEY

Example:

from FunTuple import FunTuple_Particles as Funtuple
import FunTuple.functorcollections as FC
from DaVinciMCTools import MCTruthAndBkgCat

# Load particles onto TES location for a given event cycle
input_data = get_particles("/Event/HLT2/<linename>/Particles")

# Define an algorithm that builds one-to-one relation table b/w Reco Particle -> Truth MC Particle.
# Here "options" is the DV options object (see DaVinciExamples or DaVinciTutorials)
mctruth_alg = MCTruthAndBkgCat(input_data)

variables = FC.MCPromptDecay(mctruth_alg)

tuple = Funtuple(name="MyTuple",
                 fields=...,
                 variables=variables,
                 inputs=input_data)
MCHierarchy(*, mctruth_alg: DaVinciMCTools.MCTruthAndBkgCat | None = None) FunTuple.FunctorCollection.FunctorCollection[source]

Candidate-level collection of functors yielding the hierarchy/relationships among MC particles. If no argument is provided, i.e. mctruth_alg==None, this collection creates functors that apply to MCParticles; otherwise, the functors created are assumed to apply to Particles, which will find the associated MCParticle automatically using the input truth matching.

Parameters:

mctruth_alg (DaVinciMCTools.MCTruthAndBkgCat, optional) – An instance of the MCTruthAndBkgCat helper class that is used for retrieving MC truth information of a reconstructed particle. Defaults to None.

Output variables:
  • TRUEID

  • MC_MOTHER_ID

  • MC_MOTHER_KEY

  • MC_GD_MOTHER_ID

  • MC_GD_MOTHER_KEY

  • MC_GD_GD_MOTHER_ID

  • MC_GD_GD_MOTHER_KEY

Example:

from FunTuple import FunTuple_Particles as Funtuple
import FunTuple.functorcollections as FC
from DaVinciMCTools import MCTruthAndBkgCat

# Load particles onto TES location for a given event cycle
input_data = get_particles("/Event/HLT2/<linename>/Particles")

# Define an algorithm that builds one-to-one relation table b/w Reco Particle -> Truth MC Particle.
# Here "options" is the DaVinci options object (see DaVinciExamples or DaVinciTutorials package)
mctruth_alg = MCTruthAndBkgCat(input_data)

variables = FC.MCHierarchy(mctruth_alg)

tuple = Funtuple(name="MyTuple",
                 fields=...,
                 variables=variables,
                 inputs=input_data)
MCVertexInfo(*, mctruth_alg: DaVinciMCTools.MCTruthAndBkgCat | None = None) FunTuple.FunctorCollection.FunctorCollection[source]

Candidate-level collection of functors on vertices of the related MC particle.

Parameters:
  • mctruth_alg (DaVinciMCTools.MCTruthAndBkgCat, optional) – An instance of the MCTruthAndBkgCat helper class that is used for retrieving monte-carlo truth information of a reconstructed particle. Defaults to None.

  • provided (If no argument is)

  • mctruth_alg==None (i.e.)

  • MCParticles; (this collection creates functors that apply to)

  • otherwise

  • Particles (the functors created are assumed to apply to)

:param : :param which will find the associated MCParticle automatically using the input truth matching.:

Output variables:
  • TRUEORIGINVERTEX_[X,Y,Z]

  • TRUEENDVERTEX_[X,Y,Z]

Example:

from FunTuple import FunTuple_Particles as Funtuple
import FunTuple.functorcollections as FC
from DaVinciMCTools import MCTruthAndBkgCat

# Load particles onto TES location for a given event cycle
input_data = get_particles("/Event/HLT2/<linename>/Particles")

# Define an algorithm that builds one-to-one relation table b/w Reco Particle -> Truth MC Particle.
# Here "options" is the DV options object (see DaVinciExamples or DaVinciTutorials)
mctruth_alg = MCTruthAndBkgCat(input_data)

variables = FC.MCVertexInfo(mctruth_alg)

tuple = Funtuple(name="MyTuple",
                 fields=...,
                 variables=variables,
...
)
MCPrimaries(*, mc_header: PyConf.dataflow.DataHandle) FunTuple.FunctorCollection.FunctorCollection[source]

Event-level collection of functors to store the properties (X,Y,Z,T) of all MC primary vertices.

Parameters:

mc_header (DataHandle) – The TES location of the MCHeader object, obtainded via the PyConf.reading.get_mc_header function.

Output variables:
  • MCPVX[MCPV_IDX] : F.MAP(F.VALUE_OR(F.NaN) @ F.X_COORDINATE @ F.POSITION)

  • MCPVY[MCPV_IDX] : F.MAP(F.VALUE_OR(F.NaN) @ F.Y_COORDINATE @ F.POSITION)

  • MCPVZ[MCPV_IDX] : F.MAP(F.VALUE_OR(F.NaN) @ F.Z_COORDINATE @ F.POSITION)

  • MCPVT[MCPV_IDX] : F.MAP(F.VALUE_OR(-1) @ F.MC_VTX_TIME)

  • MCPV_SIZE : F.SIZE_OF @ F.MC_ALLPVS(mc_header)

Example:

from FunTuple import FunTuple_Particles as Funtuple
import FunTuple.functorcollections as FC

mc_header = get_mc_header()

evt_variables = FC.MCPrimaries(mc_header)

tuple = Funtuple(name="MyTuple",
                 fields=...,
                 variables=...,
                 event_variables=evt_variables,
...
)
MCReconstructible(*, mcreconstructible_alg: DaVinciMCTools.MCReconstructible, extra_info: bool = False, mctruth_alg: DaVinciMCTools.MCTruthAndBkgCat | None = None) FunTuple.FunctorCollection.FunctorCollection[source]

Candidate-level collection of functors to store reconstrutible information for the MC particle.

A stable charged or neutral MCParticle can be reconstructible in one of the following categories:

-1 = No MC classification possible 0 = Outside detector acceptance 1 = In acceptance but not reconstructible 2 = Reconstructible as a Long charged track 3 = Reconstructible as a Downstream charged track 4 = Reconstructible as an Upstream charged track 5 = Reconstructible as a T charged track 6 = Reconstructible as a VELO charged track 50 = Reconstructible as a Neutral particle

Note that the support for stable neutral MCParticle is not present yet (requires separate functors). An MC decay to be reconstructible implies that all final state MCParticles are reconstructible. For more info on definition (Run1/2) see: https://lhcb-comp.web.cern.ch/analysis/davinci/v8/recrecdefinition.htm

Parameters:
  • mcreconstructible_alg (DaVinciMCTools.MCReconstructible) – An instance of the MCReconstructible helper class that builds the necessary functors.

  • extra_info (bool, optional) – If False, returns a functor to ntuple reconstructible category for charged MCParticle. If True, returns extra information e.g. hasUT, hasVelo, etc. Defaults to False.

  • mctruth_alg (DaVinciMCTools.MCTruthAndBkgCat, optional) – An instance of the MCTruthAndBkgCat helper class that is used for retrieving monte-carlo truth information of a reconstructed particle. Defaults to None.

Output variables:

  • MC_RECONSTRUCTIBLE

If extra_info:

  • MC_HASVELO

  • MC_HASVELO_AND_T

  • MC_HAST[1,2,3][ ,X,S]

  • MC_HASUT[1,2]

  • MC_ACCVELO

  • MC_ACCVELO_AND_T

  • MC_ACCT[1,2,3][ ,X,S]

  • MC_ACCUT[ ,1,2]

Example:

from FunTuple import FunTuple_MCParticles as Funtuple
import FunTuple.functorcollections as FC
from DaVinciMCTools import MCReconstructible as MCRectible

# Load MCParticles onto TES location for a given event cycle
input_data = get_particles("/Event/HLT2/MC/Particles")

# Define the helper class and get functor collection
reconstructible_alg = MCRectible()
variables = FC.MCReconstructible(reconstructible_alg)

tuple = Funtuple(name="MyTuple",
                 fields=...,
                 variables=variables,
...
)
MCReconstructed(*, mcreconstructed_alg: DaVinciMCTools.MCReconstructed, extra_info: bool = False) FunTuple.FunctorCollection.FunctorCollection[source]

Candidate-level collection of functors to store reconstrutible information for the MC particle.

A stable charged or neutral MCParticle can be reconstructed in one of the following categories:
-1 No MC classification possible (e.g. NO MC)

0 Not reconstructed 1 Reconstructed as a Long charged track 2 Reconstructed as a Downstream charged track 3 Reconstructed as an Upstream charged track 4 Reconstructed as a T charged track 5 Reconstructed as a VELO charged track 6 Reconstructed as a 2D Velo track 7 Reconstructed as a Muon track

50 Reconstructed as a Neutral particle 51 Reconstructed as a merged pi0 52 Reconstructed as a type that doesn’t fit other categories

An MC decay to be reconstructed implies that all final state MCParticles are reconstructed. For more info on definition (Run1/2) see: https://lhcb-comp.web.cern.ch/analysis/davinci/v8/recrecdefinition.htm

(since the method tuple->farray() only support vector<double> currently, all array will convert to array of floats)

Parameters:
  • mcreconstructed_alg (DaVinciMCTools.MCReconstructed) – An instance of the MCReconstructed helper class that first maps MCParticle to reconstructed track and builds the necessary functors.

  • extra_info (bool, optional) – If False, returns a functor to ntuple reconstructed category for MCParticle. If True, returns extra information e.g. TrackType, hasUT, hasVelo, TrackChi2, GhostProb, etc. Defaults to False.

Example:

from FunTuple import FunTuple_MCParticles as Funtuple
import FunTuple.functorcollections as FC
from DaVinciMCTools import MCReconstructed as MCRected

#load particles onto TES location for a given event cycle
input_data = get_particles("/Event/HLT2/MC/Particles")

#map MCParticle -> reconstructed track
mcreconstructed_alg = MCRected(input_data)
#get the variables
variables = FC.MCReconstructed(mcreconstructed_alg)

tuple = Funtuple(name="MyTuple",
                 fields=...,
                 variables=variables,
...
)
ParticleIsolation(*, isolation_alg: IsolationTools.VertexAndConeIsolation, array_indx_name: str, sumcone_invalid_value=F.NaN) FunTuple.FunctorCollection.FunctorCollection[source]

Candidate-level collection of functors on track or neutral isolation, using information from relations between sets of particles determined by a given criterion. It accesses to VertexAndConeIsolation class information.

Note::

Variables will have the following naming convention: HEAD_

Parameters:
  • isolation_alg (VertexAndConeIsolation) – An instance of the VertexAndConeIsolation helper class that is used to retrieve isolation variables.

  • array_indx_name (string) – string that helps in naming the default indx variable. This appears due to the functor arrays _DPHI, _DETA.

  • sumcone_invalid_value (optional) – Value that allows user flexibility in defining the invalid values related to the WeightedRelTable being empty. Defaults to F.NaN.

Output variables:
  • _CMULT

  • _CP[,T,X,Y,Z]

  • _P[,T,X,Y,Z]ASY

  • _DETA

  • _DPHI

Example:

import Functors as F
from FunTuple import FunTuple_Particles as Funtuple
from IsolationTools import VertexAndConeIsolation
from FunTuple.functorcollections import ParticleIsolation
from PyConf.reading import get_particles

branches = {"B": "[B0 -> ([J/psi(1S) -> tau+ mu-]CC)(K*(892)0 -> K+ pi-)]CC"}

b2ksttaumu_data = get_particles(f"/Event/HLT2/{line}/Particles")
b_cciso_data = get_particles(f"/Event/HLT2/{line}/B_{tes_long_track_iso}/Particles")

cone_isolation_b = VertexAndConeIsolation(
     name = 'TrackIsolation',
     reference_particles = b2ksttaumu_data,
     related_particles = b_cciso_data,
     cut = (F.DR<1.))

variables = FC.ParticleIsolation(isolation_alg=cone_isolation_b,
                                 array_indx_name='indx')

tuple = Funtuple(name="MyTuple",
   fields=branches,
   variables=variables,
   ...)
ConeIsolation(*, charged_cone_isolation_alg: IsolationTools.VertexAndConeIsolation, neutral_cone_isolation_alg: IsolationTools.VertexAndConeIsolation, array_indx_name: str, sumcone_invalid_value=F.NaN) FunTuple.FunctorCollection.FunctorCollection[source]

Candidate-level collection of functors on charged and neutral isolation, using information from relations between sets of particles determined by a given criterion. It accesses to VertexAndConeIsolation class information.

Note

Variables will have the following naming convention: HEAD_CC_ for charged isolation HEAD_NC_ for neutral isolation

Parameters:
  • charge_cone_isolation_alg (VertexAndConeIsolation) – An instance of the VertexAndConeIsolation helper class that is used to retrieve isolation variables for charge cone.

  • neutral_cone_isolation_alg (VertexAndConeIsolation) – An instance of the VertexAndConeIsolation helper class that is used to retrieve isolation variables for neutral cone.

  • array_indx_name (string) – string that helps in naming the default indx variable. This appears due to the functor arrays _DPHI, _DETA.

  • sumcone_invalid_value (optional) – Value that allows user flexibility in defining the invalid values related to the WeightedRelTable being empty. Defaults to F.NaN.

Output variables:

These are for both the neutral (NC) and charged (CC) cones:

  • _Max_P : charged_cone_isolation_alg.MAXCP

  • _Max_PT : charged_cone_isolation_alg.MAXCPT

  • _Min_P : charged_cone_isolation_alg.MINCP

  • _Min_PT : charged_cone_isolation_alg.MINCPT

Example:

import Functors as F
from FunTuple import FunTuple_Particles as Funtuple
from FunTuple.functorcollections import ConeIsolation
from IsolationTools import VertexAndConeIsolation
from PyConf.reading import get_particles
from GaudiKernel.SystemOfUnits import GeV

branches = {"B": "[B0 -> ([J/psi(1S) -> tau+ mu-]CC)(K*(892)0 -> K+ pi-)]CC"}

b2ksttaumu_data = get_particles(f"/Event/HLT2/{line}/Particles")
b_cciso_data = get_particles(f"/Event/HLT2/{line}/B_{tes_long_track_iso}/Particles")
b_nciso_data = get_particles(f"/Event/HLT2/{line}/B_{tes_neutrals_iso}/Particles")

charged_cone_isolation_b = VertexAndConeIsolation(
     name = 'ConeIsolation05',
     reference_particles = b2ksttaumu_data,
     related_particles = b_cciso_data,
     cut = (F.SQRT@F.DR2<0.5))

neutral_cone_isolation_b = VertexAndConeIsolation(
     name = 'ConeIsolation10',
     reference_particles = b2ksttaumu_data,
     related_particles = b_nciso_data,
     cut = (F.SQRT@F.DR2<1.))


variables = FC.ConeIsolation(charged_cone_isolation_alg=charged_cone_isolation_b,
                             neutral_cone_isolation_alg=neutral_cone_isolation_b,
                             array_indx_name='indx')

tuple = Funtuple(name="MyTuple",
   fields=branches,
   variables=variables,
   ...)
VertexIsolation(*, isolation_alg: IsolationTools.VertexAndConeIsolation) FunTuple.FunctorCollection.FunctorCollection[source]

Candidate-level collection of functors on vertex isolation, using information from relations between sets of particles determined by a given criterion. It accesses to VertexAndConeIsolation class information.

Note

Variables will have the following naming convention: VTXISO_

Parameters:
  • isolation_alg (VertexAndConeIsolation) – An instance of the VertexAndConeIsolation helper class that is used to retrieve isolation variables.

  • name (string, optional) – string that helps in distinguishing variables according to some criterion (for example if the isolation criterion are applied to charged or neutral particles). Defaults to “”.

Output variables:
  • _NParts : isolation_alg.MULT - Number of particles in CHI2 window

  • _Smallest_DELTACHI2 : isolation_alg.Smallest_DELTACHI2 - smallest delta chi2 when adding one track

  • _Smallest_CHI2 : isolation_alg.Smallest_CHI2 - smallest chi2 when adding one track

  • _Smallest_DELTACHI2_MASS : isolation_alg.Smallest_Mass_DELTACHI2 - mass of the candidate with the smallest delta chi2 when adding one track

Example:

import Functors as F
from FunTuple import FunTuple_Particles as Funtuple
from FunTuple.functorcollections import VertexIsolation
from IsolationTools import VertexAndConeIsolation
from PyConf.reading import get_particles

branches = {"B": "[B0 -> ([J/psi(1S) -> tau+ mu-]CC)(K*(892)0 -> K+ pi-)]CC"}

b2ksttaumu_data = get_particles(f"/Event/HLT2/{line}/Particles")
b_extra_one_track_cand = get_particles(f"/Event/HLT2/{line}/B_{tes_long_one_track_iso}/Particles")
b_extra_two_tracks_cand = get_particles(f"/Event/HLT2/{line}/B_{tes_long_two_tracks_iso}/Particles")

vertex_isolation_b_one_track = VertexAndConeIsolation(
     name = 'OneTrack',
     reference_particles = b2ksttaumu_data,
     related_particles = b_extra_one_track_cand,
     cut = (F.CHI2@F.F.FORWARDARG1<9.))

vertex_isolation_b_two_tracks = VertexAndConeIsolation(
     name = 'TwoTracks',
     reference_particles = b2ksttaumu_data,
     related_particles = b_extra_two_tracks_cand,
     cut = (F.CHI2@F.F.FORWARDARG1<15.))

variables_ot = FC.VertexVariables(isolation_alg=vertex_isolation_b_one_track)
variables_tt = FC.VertexVariables(isolation_alg=vertex_isolation_b_two_tracks)

tuple = Funtuple(name="MyTuple",
   fields=branches,
   variables=variables_ot+variables_tt,
   ...)
NeutralCaloInfo(extra_info: bool = False) FunTuple.FunctorCollection.FunctorCollection[source]

Candidate-level collection of functors on neutral calorimeter hypotheses information. Useful for example for neutral particle reconstruction efficiency or neutral PID studies.

Parameters:

extra_info (bool, optional) – If True, adds Neutral CellIDs bitwise information, row and column.

Output variables:
  • CaloTrackMatchChi2 : 2D chi2 for Track-CaloCluster matching

  • CaloNeutralShowerShape : 2nd order moment of the cluster

  • CaloClusterMass : CALO cluster Mass obtained from MergedPi0Alg

  • CaloNeutralEcalEnergy : The ECAL cluster energy associated to the neutral CaloHypo

  • CaloNeutral1To9EnergyRatio : Fraction of highest 1x1 cell to 3x3 cells forming the cluster for the Neutral CaloHypo

  • CaloNeutral4To9EnergyRatio : Fraction of highest 2x2 cell to 3x3 cells forming the cluster for the Neutral CaloHypo

  • CaloNeutralHcal2EcalEnergyRatio : The Hcal/Ecal energy ratio associated to the neutral calorimeter hypotheses

  • CaloNumSaturatedCells : Number of saturated cells

  • CaloNeutralArea : Area code of Neutral CellID

If extra_info is set to True, it adds also:
  • CaloNeutralID : Bitwise information (32 bits) of all Neutral seed CellIDs

  • CaloNeutraRow : Row code of Neutral CellID

  • CaloNeutraCol : Column code of Neutral CellID

ChargedCaloInfo(extra_info: bool = False) FunTuple.FunctorCollection.FunctorCollection[source]

Candidate-level collection of functors on charged calorimeter hypotheses information. Useful for example for charged particle reconstruction efficiency or charged PID studies.

Parameters:

extra_info (bool, optional) – If True, adds extra information, including CaloCluster and CaloHypo (photon and electron) CellIDs related ones.

Output variables:
  • INECAL : In Ecal acceptance

  • INHCAL : In Hcal acceptance

  • INBREM : In Brem acceptance: the track extrapolation from upstream of the magnet is in Ecal acceptance.

  • HASBREM : Has non-zero brem momentum-recovery energy available

  • HASBREMADDED : Has non-zero brem momentum-recovery energy added by BremAdder

  • BREMPIDE : Brem-based DLL for electron-ID

  • ECALPIDE : Ecal-based DLL for electron-ID

  • ECALPIDMU : Ecal-based DLL for muon-ID

  • HCALPIDE : Hcal-based DLL for electron-ID

  • HCALPIDMU : Hcal-based DLL for muon-ID

  • CLUSTERAREA : Area code of CaloCluster CellID

  • ELECTRONAREA : Area code of CaloHypo (electron) CellID

  • BREMHYPOAREA : Area code of CaloHypo (photon) CellID

If extra_info is set to True, it adds also:
  • CLUSTERID : All significant bits representation of CellID (32bits), i.e. CellID.all(), associated to CaloCluster of track match

  • CLUSTERROW : Row code of CaloCluster CellID

  • CLUSTERCOL : Column code of CaloCluster CellID

  • CLUSTERMATCH_CHI2 : 2D chi2 for Track/CaloCluster matching (neutral/charged)

  • ELECTRONSHOWEREOP : Ratio of electron energy/momentum with track-based cell selection

  • ELECTRONSHOWERDLL : Summed per-cell E/p DLL (electron versus pion) with track-based cell selection and energy estimation

  • ELECTRONMATCH_CHI2 : 3D chi2 for Track/CaloHypo(e) matching (charged)

  • ELECTRONENERGY : Cluster energy associated to CaloHypo (charged)

  • ELECTRONID : All significant bits representation of CellID (32bits), i.e. CellID.all(), associated to CaloHypo seed (electron hypo).

  • ELECTRONROW : Row code of CaloHypo (electron) CellID

  • ELECTRONCOL : Column code of CaloHypo (electron) CellID

  • BREMENERGY : Brem momentum-recovery energy

  • BREMBENDCORR : Correction factor accounting for bending biases in track due to brem

  • BREMHYPOID : All significant bits representation of CellID (32bits), i.e. CellID.all(), for CaloHypo (photon) associated to track for brem recovery

  • BREMHYPOROW : Row code of CaloHypo (photon) CellID

  • BREMHYPOCOL : Column code of CaloHypo (photon) CellID

  • BREMHYPOMATCH_CHI2 : 2D chi2 of CaloHypo (photon) associated to track for brem recovery

  • BREMHYPOENERGY : Energy of CaloHypo (photon) associated to track for brem recovery

  • BREMHYPODELTAX : Test statistic of being first-state like of CaloHypo (photon) for brem recovery

  • BREMTRACKBASEDENERGY : Track-based brem energy determination

  • HCALEOP : Hcal energy deposit over momentum (track)

ParticleID(extra_info: bool = False) FunTuple.FunctorCollection.FunctorCollection[source]

Candidate-level collection of functors on global particle identification variables.

Parameters:

extra_info (bool, optional) – If True, store DLLs in addition. Defaults to False.

Output variables:
  • PARTICLE_ID : F.PARTICLE_ID

  • PROBNN_[D,E,GHOST,K,MU,P,PI] : F.PROBNN_[D,E,GHOST,K,MU,P,PI]

RecSummary(detector_info: bool = False) FunTuple.FunctorCollection.FunctorCollection[source]

Event-level collection providing access to the RecSummary record. Provides event multiplicities. The list of quantities is at https://gitlab.cern.ch/lhcb/LHCb/-/blob/master/Event/RecEvent/include/Event/RecSummary.h.

Parameters:

detector_info (bool, optional) – If True, stores hit multiplicities in subdetectors. Defaults to False.

Output variables:
  • nPVs

  • n[Long,Downstream,Upstream,Velo,T,Back]Tracks

  • [e,h]CalTot

  • nRich[1,2]Hits

  • n[VP,UT,FT,Ecal]Clusters

Example:

import FunTuple.functorcollections as FC
from DaVinci import options

variables = FC.RecSummary()
# ...

tuple = FuntupleEvent(name="MyTuple",
                      fields=...,
                      variables=variables,
                      ...)
FlavourTaggingResults(all_taggers: PyConf.dataflow.DataHandle, tagger_name_list: list = ['Run2_SSPion', 'Run2_SSKaon', 'Run2_SSProton', 'Run2_OSKaon', 'Run2_OSElectron', 'Run2_OSMuon', 'Run2_OSVertexCharge'], prefix: str = '') FunTuple.FunctorCollection.FunctorCollection[source]

Candidate-level collection of functors that access FlavourTagging results.

Parameters:
  • all_taggers – a FlavourTags object that holds the information of all relevant taggers

  • tagger_name_list (list, optional) – a list of all the names (as specified in LHCb::Tagger) of the taggers that are requested as branches. Default is all Run2 taggers.

  • prefix (str, optional) – additional prefix to output functors. Defaults to “”.

Output variables:
  • _Dec : F.TaggingDecision for each tagger specified

  • _Omega : F.TaggingMistag for each tagger specified

  • _MVA : F.TaggingMVAOutput for each tagger specified

Example:

from FunTuple import FunTuple_Particles as Funtuple
import FunTuple.functorcollections as FC
from Hlt2Conf.flavourTagging import run2_all_taggers

line_data = get_particles(f"/Event/HLT2/{line_name}/Particles")
all_tagging = run2_all_taggers(line_data)

variables = {"B": FC.FlavourTaggingResults(all_tagging) }

tuple = Funtuple(
    name="Tuple",
    tuple_name="DecayTree",
    fields=fields,
    variables=variables,
    inputs=line_data,
)