Belle II Software development
Chi2MCTrackMatcherModule.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#include <tracking/modules/mcMatcher/Chi2MCTrackMatcherModule.h>
10
11// datastore types
12#include <framework/datastore/StoreArray.h>
13#include <mdst/dataobjects/Track.h>
14#include <tracking/dataobjects/RecoTrack.h>
15#include <mdst/dataobjects/MCParticle.h>
16#include <framework/gearbox/Const.h>
17#include <Eigen/Dense>
18#include <iostream>
19
20// ROOT
21#include <TVectorD.h>
22#include <TMatrixD.h>
23
24
25using namespace Belle2;
26
27REG_MODULE(Chi2MCTrackMatcher);
28
29
30
31// Constructor
33{
34 // Set module properties
35 setDescription(R"DOC(Monte Carlo matcher using the helix parameters for matching by chi2-method)DOC");
36 // setPropertyFlags(c_ParallelProcessingCertified);
37 // set module parameters
38 addParam("CutOffs",
40 "Defines the chi2-cut-off values for each charged particle. The cut-offs define the maximal size for a matching pair candidate`s chi2 value. If a matching pair candidate has a chi2 value smaller than the cut-off a relation is set. Defaults determined from small study on events with trivial matching. The pdg order is [11,13,211,2212,321,1000010020].",
42//,defaultCutOffs);
43 addParam("linalg",
45 "Parameter to switch between ROOT and Eigen, to invert the covariance5 matrix. ROOT has shown a shorter runtime therefore its recommended. false: ROOT is used; true: Eigen is used",
46 false);
47 addParam("MCRecoTracksArrayName",
49 "Name of the StoreArray of MC RecoTracks which will be related to the PR RecoTracks, "
50 "using the relations to Belle2::MCParticles and Belle2::Tracks",
52 addParam("PRRecoTracksArrayName",
54 "Name of the StoreArray of PR RecoTracks which will be related to the MC RecoTracks, "
55 "using the relations to Belle2::MCParticles and Belle2::Tracks",
57}
58
60{
61 // Check if there are MC Particles
63 m_Tracks.isRequired();
64 m_Tracks.registerRelationTo(m_MCParticles);
65
67 mcRecoTracks.isOptional();
69 prRecoTracks.isOptional();
70
71 // if both arrays exist we want a relation between them
72 if (mcRecoTracks.isValid() and prRecoTracks.isValid()) {
73 prRecoTracks.registerRelationTo(mcRecoTracks);
74 mcRecoTracks.registerRelationTo(prRecoTracks);
77 }
78
79}
80
82{
83 m_event += 1;
84 B2DEBUG(27, "m_event = " << m_event);
85 // suppress the ROOT "matrix is singular" error message
86 auto default_gErrorIgnoreLevel = gErrorIgnoreLevel;
87 gErrorIgnoreLevel = 5000;
88 // get StoreArray length
89 int nTracks = m_Tracks.getEntries();
90 int nMCParticles = m_MCParticles.getEntries();
91 // tracks number of tracks for Debug statistics
92 m_trackCount += nTracks;
93 // check if there are m_Tracks and m_MCParticles to match to
94 if (not nMCParticles or not nTracks) {
95 // Cannot perform matching
96 return;
97 }
98 // compare all tracks with all mcParticles in event
99 for (auto& track : m_Tracks) {
100 // test for existing relations
101 if (track.getRelated<MCParticle>()) {
102 B2DEBUG(27, "Relation already set continue with next track");
103 continue;
104 }
105 // initialize minimal chi2 variables
106 MCParticle* mcPart_matched = nullptr;
107 double chi2Min = std::numeric_limits<double>::infinity();
108
109 // loop over m_MCParticles and calculate Chi2 for each track mcparticle pair, to fine minimal chi2
110 for (auto& mcParticle : m_MCParticles) {
111 // Check if current particle is a charged stable particle
112 const Const::ParticleType mcParticleType(std::abs(mcParticle.getPDG()));
113 if (not Const::chargedStableSet.contains(mcParticleType)) {
114 continue;
115 }
116 // get trackfitresult of current track with closest mass to current mcparticle Type
117 auto trackFitResult = track.getTrackFitResultWithClosestMass(mcParticleType);
118 TMatrixD Covariance5 = trackFitResult->getCovariance5();
119 // statistic variable counting number of covariance matrices
121 // check if matrix is invertable
122 double det = Covariance5.Determinant();
123 if (det == 0.0) {
124 // statistic variable counting number of not invertable covariance matrices
126 //Covariance5.Print();
127 continue;
128 }
129 // generate helix for current mc particle
130 double charge_sign = 1;
131 if (mcParticle.getCharge() < 0) { charge_sign = -1;}
132 UncertainHelix helix = trackFitResult->getUncertainHelix();
133 auto b_field = BFieldManager::getField(helix.getPerigee()).Z() / Belle2::Unit::T;
134 auto mcParticleHelix = Helix(mcParticle.getVertex(), mcParticle.getMomentum(), charge_sign, b_field);
135 // initialize the current chi2
136 double chi2Cur;
137 // Check which linear algebra system should be used and calculate chi2Cur
138 if (not m_param_linalg) {
139 // Eigen
140 // build difference vector
141 Eigen::VectorXd delta(5);
142 delta(0) = trackFitResult->getD0() - mcParticleHelix.getD0();
143 delta(1) = trackFitResult->getPhi0() - mcParticleHelix.getPhi0();
144 delta(2) = trackFitResult->getOmega() - mcParticleHelix.getOmega();
145 delta(3) = trackFitResult->getZ0() - mcParticleHelix.getZ0();
146 delta(4) = trackFitResult->getTanLambda() - mcParticleHelix.getTanLambda();
147 // build convariance5 Eigen matrix
148 Eigen::MatrixXd covariance5Eigen(5, 5);
149 for (int i = 0; i < 5; i++) {
150 for (int j = 0; j < 5; j++) {
151 covariance5Eigen(i, j) = Covariance5[i][j];
152 }
153 }
154 //calculate chi2Cur
155 chi2Cur = ((delta.transpose()) * (covariance5Eigen.inverse() * delta))(0, 0);
156 } else {
157 // ROOT
158 // calculate difference vector
159 TMatrixD delta(5, 1);
160 delta[0][0] = trackFitResult->getD0() - mcParticleHelix.getD0();
161 delta[1][0] = trackFitResult->getPhi0() - mcParticleHelix.getPhi0();
162 delta[2][0] = trackFitResult->getOmega() - mcParticleHelix.getOmega();
163 delta[3][0] = trackFitResult->getZ0() - mcParticleHelix.getZ0();
164 delta[4][0] = trackFitResult->getTanLambda() - mcParticleHelix.getTanLambda();
165 // invert Covariance5 matrix
166 Covariance5.InvertFast();
167 // transpose difference vector
168 TMatrixD deltaT = delta;
169 deltaT.T();
170 // calculate chi2Cur
171 chi2Cur = ((deltaT) * (Covariance5 * delta))[0][0];
172 }
173 // check if chi2Cur is smaller than the so far found minimal chi2Min
174 if (chi2Min > chi2Cur) {
175 chi2Min = chi2Cur;
176 mcPart_matched = &mcParticle;
177 }
178 }
179 // check if any matching candidate was found
180 if (not mcPart_matched) {
182 m_noMatchCount += 1;
183 continue;
184 }
185 // initialize Cut Off
186 double cutOff = 0;
187 // fill cut off with value
188 // find PDG for mcParticle with minimal chi2, because this determines individual cutOff
189 int mcMinPDG = abs(mcPart_matched->getPDG());
190 if (mcMinPDG == Belle2::Const::electron.getPDGCode()) {
191 cutOff = m_param_CutOffs[0];
192 } else if (mcMinPDG == Const::muon.getPDGCode()) {
193 cutOff = m_param_CutOffs[1];
194 } else if (mcMinPDG == Const::pion.getPDGCode()) {
195 cutOff = m_param_CutOffs[2];
196 } else if (mcMinPDG == Const::proton.getPDGCode()) {
197 cutOff = m_param_CutOffs[3];
198 } else if (mcMinPDG == Const::kaon.getPDGCode()) {
199 cutOff = m_param_CutOffs[4];
200 } else if (mcMinPDG == Const::deuteron.getPDGCode()) {
201 cutOff = m_param_CutOffs[5];
202 } else {
203 B2WARNING("The pdg for minimal chi2 was not in chargedstable particle list: MinPDG = " << mcMinPDG);
204 continue;
205 }
206 if (chi2Min < cutOff) {
207 track.addRelationTo(mcPart_matched);
208 // set relations between MC and PR RecoTracks, needed by tracking validation
209 RecoTrack* mcRecoTrack = mcPart_matched->getRelated<RecoTrack>(m_MCRecoTracksArrayName);
210 RecoTrack* prRecoTrack = track.getRelated<RecoTrack>(m_PRRecoTracksArrayName);
211 if (mcRecoTrack and prRecoTrack) {
212 prRecoTrack->setMatchingStatus(RecoTrack::MatchingStatus::c_matched);
213 prRecoTrack->addRelationTo(mcRecoTrack);
214 mcRecoTrack->addRelationTo(prRecoTrack);
215 prRecoTrack->addRelationTo(mcPart_matched);
216 mcPart_matched->addRelationTo(prRecoTrack);
217 }
218 } else {
219 m_noMatchCount += 1;
220 }
221 }
222 // reset ROOT Error Level to default
223 gErrorIgnoreLevel = default_gErrorIgnoreLevel;
224}
226{
227 B2DEBUG(21, "For " << m_noMatchCount << "/" << m_trackCount << "(" << m_noMatchCount * 100 / m_trackCount <<
228 "%) tracks no relation was found.");
229 B2DEBUG(21, "For " << m_noMatchingCandidateCount << "/" << m_trackCount << "(" << m_noMatchingCandidateCount * 100 / m_trackCount <<
230 "%) tracks got no matching candidate");
231 B2DEBUG(21, "" << m_notInvertableCount << "/" << m_covarianceMatrixCount << "(" << m_notInvertableCount * 100 /
232 m_covarianceMatrixCount << "%) covariance matrices where not invertable.");
233}
234
Helix parameter class.
Definition: Helix.h:48
bool m_param_linalg
Parameter: Possibility to switch between ROOT and Eigen for inversion of the covariance matrix false:...
Chi2MCTrackMatcherModule()
Constructor: Sets the description, the properties and the parameters of the module.
std::string m_PRRecoTracksArrayName
Variable: Name of the RecoTracks StoreArray of PR tracks.
int m_noMatchCount
Variable: Counts the number of tracks without a relation.
std::vector< double > m_param_CutOffs
Parameter : Defines the chi2 cut-off values for each charged particle.
void initialize() override
Register input and output data.
void event() override
Do matching for each event.
int m_event
Variable: Counts total event number.
std::string m_MCRecoTracksArrayName
Variable: Name of the RecoTracks StoreArray of MC tracks.
void terminate() override
Provides debug statistics.
int m_noMatchingCandidateCount
Variable: Counts the total number of tracks without a possible MCParticle as matching candidate.
int m_notInvertableCount
Variables for Debug module statistics.
int m_trackCount
Variable: Counts the total number of tracks to calculate a percentage feedback.
int m_covarianceMatrixCount
Variable: Counts the total number of Covaraince5 matrices.
StoreArray< Track > m_Tracks
Variable: Makes m_Tracks available in whole class.
StoreArray< MCParticle > m_MCParticles
Variable: Makes m_MCParticles available in whole class.
The ParticleType class for identifying different particle types.
Definition: Const.h:408
static const ChargedStable muon
muon particle
Definition: Const.h:660
static const ParticleSet chargedStableSet
set of charged stable particles
Definition: Const.h:618
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
static const ChargedStable proton
proton particle
Definition: Const.h:663
static const ChargedStable kaon
charged kaon particle
Definition: Const.h:662
static const ChargedStable electron
electron particle
Definition: Const.h:659
static const ChargedStable deuteron
deuteron particle
Definition: Const.h:664
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
void setMatchingStatus(MatchingStatus matchingStatus)
Set the matching status (used by the TrackMatcher module)
Definition: RecoTrack.h:835
void addRelationTo(const RelationsInterface< BASE > *object, float weight=1.0, const std::string &namedRelation="") const
Add a relation from this object to another object (with caching).
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
bool isValid() const
Check wether the array was registered.
Definition: StoreArray.h:288
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
bool registerRelationTo(const StoreArray< TO > &toArray, DataStore::EDurability durability=DataStore::c_Event, DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut, const std::string &namedRelation="") const
Register a relation to the given StoreArray.
Definition: StoreArray.h:140
This class represents an ideal helix in perigee parameterization including the covariance matrix of t...
static const double T
[tesla]
Definition: Unit.h:120
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
static void getField(const double *pos, double *field)
return the magnetic field at a given position.
Definition: BFieldManager.h:91
Abstract base class for different kinds of events.