7.5.1. Tree Fitter#

The TreeFitter is a global fitting tool to simultaneously fit an entire decay chain. It tries to find the best decay vertex positions and momenta for the given decay hypothesis. Final state particles will be constrained to their masses, meaning no energy is extracted from the fit but pdg-mass and momentum measurement are used to determine it (same as track fitting procedure). The fit extracts the 4-momenta (->including energy) and decay vertex positions of intermediate particles. If an intermediate particle has a flight length less than \(1~\mu m\) - for example a strongly decaying hadronic resonance - its decay and production vertex will be merged. The production vertex is the decay vertex of the mother particle. Note that the other fitters do not do this.

The fit tries to find the solutions to a system of equations or at least the solution closest to the measurement given the measurement uncertainty and noise. The distance between (linearised) hypothesis and measurement is calculated using \(\chi^2 = \sum \frac{h-m}{\sigma}\) the final goodness of the fit then can be evaluated using the p-value defined as the incomplete (lower = 1 - upper) gamma function of the \(\chi^2\) and degrees of freedom \(r\) divided by two: \(1 -\int^{ \frac{\chi^2}{2} }_{0} t^{r/2-1} e^{-t} dt~/~\Gamma(\frac{r}{2}).\)

This quantity denotes the probability the observed \(\chi^2\) exceeds the expectations for a correct model by chance (see ROOT documentation TMath::Prob(Double_t chi2,Int_t ndf)).

The number of degrees of freedom is the number of equations (constraints) minus the number of parameters extracted (vertices, momenta, …) and has to be larger than zero, otherwise the fit will display a warning and do nothing. The constraints can be counted as following:

Particle

parameters

constraints

track

px, py, pz

5 (helix)

photon

px, py, pz

3 (see paper)

Intermediate

px, py, pz, E, x, y, z, tau

4 (kinematic) + 3 (geometric) + 1 (mass)

\(\tau\) is the decay length. For intermediates it depends, generally there are 4 equations from the kinematic constraint, if they have a flight length longer than \(1~\mu m\) AND a well defined production vertex (a mother particle or an origin constraint) we extract the flight length and can use a geometric constraint, which the fit will automatically do unless configured otherwise. If the particle shares the decay vertex with the mother we will only extract this vertex parameters once in the fit. A mass constraint will only be used if the user configures the fit to do so. So for example for \(J/\psi \rightarrow l^+ l^-\) the calculation is the following: \(r = 5 (track) + 5 (track) + 4 (kinematic) \\- 3 (track) - 3 (track) \\- 7 (J/\psi:~px, py, pz, E, x, y, z) = 1.\)

While for \(B \rightarrow D(K \pi \pi^0)\pi\): \(r = 3\dot 5 (track) + 2 \cdot 3 (gamma) \\+ 3 \cdot 4 (kinematic) + 3 (D:geometric) \\- 6 (B:~x,y,z,px,py,pz) - 7 (D:~x,y,z,tau,px,py,pz) \\- 3 (\pi^0:~px,py,pz) - 3 \cdot 3 (track) \\- 2 \cdot 3 (gamma) = 5.\)

Citation#

Krohn, J.-F. et al. Nucl.Instrum.Meth.A 976 (2020) 164269 https://doi.org/10.1016/j.nima.2020.164269

Usage#

The user reconstructs the desired decay chain and then passes the head of the tree to the vertex fitter. So when you are reconstructing \(B->D(K \pi)\pi\) the list name of the B-meson has to be passed and the fit fits everything down the stream. The vertex fitter has a convenience function: vertex.treeFit.

Hint

The TreeFit will most likely change the kinematic properties of the decay head (mother) and, if the option updateAllDaughters is turned on, also of the daughter particles. The variables of daughters before the TreeFit can be accessed with the meta-variable originalDaughter. For the variables of the mother particle, one can store the quantities in the extraInfo field by running variablesToExtraInfo before the TreeFit and then accessing them via the variable extraInfo.

Troubleshooting and FAQ#

  • The TreeFitter is extremely sensitive to how the initial covariance matrix is initialised, we haven’t found a good general way to do it. If you see that the signal efficiency of your fit is below 95% this might be the problem. Feel free to report this. But keep in mind garbage in -> garbage out.

  • Vertex fitting MC particles (ie. particle lists created with modularAnalysis.fillParticleListsFromMC) can lead to crashes. MC particles do not need to be and should not be vertexed. If you wish to study the vertex resolution please fit the reconstructed particles and compare the vertex variables with those of the truth matched MC values.

See also

Ask for help with TreeFitter (or ask questions about it) at questions.belle2.org. Please tag your questions with TreeFitter.

Parameters of the convenience function#

Import and use the convenience function:

vertex.treeFit(list_name, conf_level=0.001, massConstraint=[], ipConstraint=False, updateAllDaughters=False, customOriginConstraint=False, customOriginVertex=[0.001, 0, 0.0116], customOriginCovariance=[0.0048, 0, 0, 0, 0.003567, 0, 0, 0, 0.04], originDimension=3, treatAsInvisible='', ignoreFromVertexFit='', path=None)[source]

Perform a Tree Fitter fit for each Particle in the given ParticleList.

Example

An example of usage for the decay chain \(B^0\to\pi^+\pi^-\pi^0\) is the following:

reconstructDecay('pi0:A -> gamma:pi0 gamma:pi0', '0.130 < InvM < 0.14', path=mypath)
reconstructDecay('B0:treefit -> pi+:my pi-:my pi0:A ', '', path=mypath)
treeFit('B0:treefit', ipConstraint=True, path=mypath)

Note

The Particle object loaded from the KLMCluster does not have proper covariance matrix yet (by release-06 at least). The performance of TreeFit with these particles is not guaranteed. Alternatively, one can perform the TreeFit ignoring these particles from the vertex fit with the option ignoreFromVertexFit.

Parameters
  • list_name (str) – name of the input ParticleList

  • conf_level (float) – minimum value of the confidence level to accept the fit. Setting this parameter to -1 selects all particle candidates. The value of 0 rejects the particle candidates with failed fit.

  • massConstraint (list(int) or list(str)) – list of PDG ids or Names of the particles which are mass-constrained Please do not mix PDG id and particle names in massConstraint list.

  • ipConstraint (bool) – constrain head production vertex to IP (x-y-z) constraint

  • customOriginConstraint (bool) – use a custom origin vertex as the production vertex of your particle. This is useful when fitting D*/D without wanting to fit a B but constraining the process to be B-decay-like. (think of semileptonic modes and stuff with a neutrino in the B decay).

  • customOriginVertex (list(float)) – 3d vector of the vertex coordinates you want to use as custom origin. Default numbers are taken for B-mesons

  • customOriginCovariance (list(float)) – 3x3 covariance matrix for the custom vertex (type: vector). Default numbers extracted from generator distribution width of B-mesons.

  • updateAllDaughters (bool) – if true the entire tree will be updated with the fitted values for momenta and vertex position. Otherwise only the momenta of the head of the tree will be updated, however for all daughters we also update the vertex position with the fit results as this would otherwise be set to {0, 0, 0} contact us if this causes any hardship/confusion.

  • originDimension (int) – If the origin or IP constraint (customOriginVertex or ipConstraint) are used, this specifies the dimension of the constraint (3D or 2D).

  • treatAsInvisible (str) – Decay string to select one particle that will be treated as invisible in the fit.

  • ignoreFromVertexFit (str) – Decay string to select particles that will be ignored to determine the vertex position.

  • path (basf2.Path) – modules are added to this path

All other parameters are the same as for the module, see below. The user has to set a parameter confidenceLevel that corresponds to the p-value as calculated above. Fits with values below the set value will be discarded as junk.

(Advanced) module configuration#

The following are parameters as they would be passed to THE MODULE. The convenience function has slightly different names to make it more pythonic, see the doc of the convenience function. The meaning is the same. To use the module add it to the path:

list_name = 'B0:list_name'
path.add_module('TreeFitter',
    particleList=list_name
    ...
    )
# additionally set a name for the module in the path, not necessary
for m in path.modules():
  if "TreeFitter" ==  m.name():
    m.set_name(f'TreeFitter_{list_name}')
TreeFitter

Tree Fitter module. Performs simultaneous fit of all vertices in a decay chain. Can also be used to just fit a single vertex.

Package

analysis

Library

libTreeFitter.so

Required Parameters
  • particleList (str)

    Type::[string]. Input mother of the decay tree to fit. For example ‘B0:myB0particleList’.

Parameters
  • autoSetGeoConstraintAndMergeVertices (bool, default=True)

    Type::bool, shall vertices of strong resonance be merged with their mothers? Can the particles vertex be constraint geometrically?

  • confidenceLevel (float, default=0.0)

    Type::[double]. Confidence level to accept fitted decay tree. Candidates with < confidenceLevel will be removed from the particle list! Typical Values: -1: keep all particle candidates, 0: remove all that fail the fit, 0.001: standard cut, 0.1: (too) tight cut. Optimise using a figure of merit (for example S/(sqrt{S+B}) ) for your analysis.

  • convergencePrecision (float, default=0.01)

    Type::[double]. Fractional upper limit for chi2 fluctuations to accept result. Larger value = less signal rejection but also less background rejection. Optimized for FOM on different topologies - don’t touch unless you REALLY want this.

  • customOriginConstraint (bool, default=False)

    Type::[bool]. Use a custom vertex as the production point of the highest hierarchy particle (register this as the mother of the list you specify). Like the beam constraint but you can specify the position its covariance yourself.

  • customOriginCovariance (list(float), default=[0.0048, 0.0, 0.0, 0.0, 0.003567, 0.0, 0.0, 0.0, 0.04])

    Type::[double]. List vertex covariance elements used in the custom origin constraint (as a vector). Default is meant for B0 decays and is taken from 100k generated B0 to mumu events.

  • customOriginVertex (list(float), default=[0.001, 0.0, 0.0116])

    Type::[double]. List of vertex coordinates to be used in the custom origin constraint.

  • expertBeamConstraintPDG (int, default=0)

    Type int, default 0. The 4-momentum of particles with the given PDG will be constrained to the 4-momentum of the initial e+e- system.

  • expertMassConstraintType (int, default=0)

    Type::[int]. False(0): use particles parameters in mass constraint; True: use sum of daughter parameters for mass constraint. WAARNING not even guaranteed that it works.

  • expertRemoveConstraintList (list(str), default=[])

    Type::[string]. List of constraints that you do not want to be used in the fit. WARNING don’t use if you don’t know exactly what it does.

  • expertUseReferencing (bool, default=True)

    Type::[bool]. Use the Extended Kalman Filter. This implementation linearises around the previous state vector which gives smoother convergence.

  • geoConstraintList (list(int), default=[])

    Type::[int], if ‘autoSetGeoConstraintAndMergeVertices==False’ you can manually set the particles that will be geometrically constrained here.

  • ignoreFromVertexFit (str, default=’’)

    Type::[string]. Decay string to select particles that will be ignored to determine the vertex position while kept for kinematics determination.

  • inflationFactorCovZ (int, default=1)

    Inflate the covariance of the beamspot by this number so that the 3d beam constraint becomes weaker in Z.And: thisnumber->infinity : dim(beamspot constr) 3d->2d.

  • ipConstraint (bool, default=False)

    Type::[bool]. Use the IP as the origin of the tree. This registers an internal IP particle as the mother of the list you give. Or in other words forces the PRODUCTION vertex of your particle to be the IP and its covariance as specified in the database.

  • massConstraintList (list(int), default=[])

    Type::[int]. List of particles to mass constrain with int = pdg code. Note that the variables ‘M’: fit result for the particle and ‘InvM’: calculated from the daughter momenta, will look different (especially if you don’t update the daughters!).

  • massConstraintListParticlename (list(str), default=[])

    Type::[string]. List of particles to mass constrain with string = particle name.

  • originDimension (int, default=3)

    Type int, default 3. If origin or ip constraint used, specify the dimension of the constraint 3->x,y,z; 2->x,y. This also changes the dimension of the geometric constraints! So you might want to turn them off for some particles. (That means turn auto off and manually on for the ones you want to constrain)

  • sharedVertexList (list(int), default=[])

    Type::[int], if ‘autoSetGeoConstraintAndMergeVertices==False’ you can manually set the particles that share the vertex with their mother here.

  • treatAsInvisible (str, default=’’)

    Type::[string]. Decay string to select one particle that will be ignored in the fit.

  • updateAllDaughters (bool, default=False)

    Type::[bool]. Update all daughters (vertex position and momenta) in the tree. If not set only the 4-momenta for the head of the tree will be updated. We also update the vertex position of the daughters regardless of what you put here, because otherwise the default when the particle list is created is {0,0,0}.

There are more parameters which are based on optimisations we made. It is strongly recommended to not touch them.