Belle II Software development
FillTrackFitNtupleModule.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/trackingPerformanceEvaluation/FillTrackFitNtupleModule.h>
10
11#include <mdst/dataobjects/HitPatternCDC.h>
12#include <mdst/dataobjects/HitPatternVXD.h>
13
14#include <map>
15
16using namespace Belle2;
17
18//-----------------------------------------------------------------
19// Register the Module
20//-----------------------------------------------------------------
21REG_MODULE(FillTrackFitNtuple);
22
24 Module()
25{
26
27 setDescription("This module fills a ntuple with tracking variables under different hypotheses");
28
29 addParam("outputFileName", m_rootFileName, "Name of output root file.",
30 std::string("TrackingPerformanceEvaluation_output.root"));
31 addParam("TracksName", m_TracksName, "Name of Track collection.", std::string(""));
32 addParam("RecoTracksName", m_RecoTracksName, "Name of RecoTrack collection.", std::string("RecoTracks"));
33 addParam("ParticleHypothesis", m_ParticleHypothesis, "Particle Hypothesis used in the track fit.", int(211));
34
35}
36
38{
39 // Tracks, RecoTracks needed for this module
41 m_Tracks.isRequired();
42
43 //set the ROOT File
44 m_rootFilePtr = new TFile(m_rootFileName.c_str(), "RECREATE");
45
46 //now create ntuples
47 TString var_list("evt:run:exp:prod:nhits:ncdc:npxd:nsvd:seed_x:");
48 var_list += ("seed_y:seed_z:seed_px:seed_py:seed_pz:seed_p:seed_pt:seed_theta:seed_phi:seed_charge:");
49 var_list += ("nhits_pi:ncdc_pi:npxd_pi:nsvd_pi:nhits_k:ncdc_k:npxd_k:nsvd_k:nhits_p:ncdc_p:npxd_p:nsvd_p:nhits_d:ncdc_d:npxd_d:nsvd_d:");
50 var_list += ("flag_pi:flag_k:flag_p:flag_d:");
51 var_list += ("trk_x_pi:trk_y_pi:trk_z_pi:trk_px_pi:trk_py_pi:trk_pz_pi:trk_p_pi:trk_pt_pi:trk_theta_pi:trk_phi_pi:");
52 var_list += ("trk_charge_pi:trk_chi2_pi:trk_ndf_pi:trk_pvalue_pi:nfailed_pi:");
53 var_list += ("trk_x_k:trk_y_k:trk_z_k:trk_px_k:trk_py_k:trk_pz_k:trk_p_k:trk_pt_k:trk_theta_k:trk_phi_k:");
54 var_list += ("trk_charge_k:trk_chi2_k:trk_ndf_k:trk_pvalue_k:nfailed_k:");
55 var_list += ("trk_x_p:trk_y_p:trk_z_p:trk_px_p:trk_py_p:trk_pz_p:trk_p_p:trk_pt_p:trk_theta_p:trk_phi_p:");
56 var_list += ("trk_charge_p:trk_chi2_p:trk_ndf_p:trk_pvalue_p:nfailed_p:");
57 var_list += ("trk_x_d:trk_y_d:trk_z_d:trk_px_d:trk_py_d:trk_pz_d:trk_p_d:trk_pt_d:trk_theta_d:trk_phi_d:");
58 var_list += ("trk_charge_d:trk_chi2_d:trk_ndf_d:trk_pvalue_d:nfailed_d:");
59 var_list += ("cdcf_pi:cdcl_pi:svdf_pi:svdl_pi:");
60 var_list += ("cdcf_k:cdcl_k:svdf_k:svdl_k:");
61 var_list += ("cdcf_p:cdcl_p:svdf_p:svdl_p:");
62 var_list += ("cdcf_d:cdcl_d:svdf_d:svdl_d");
63 m_n_MultiParticle = new TNtuple("nMultiParticle", "ntuple for multi hyp particle", var_list);
64}
65
67{
68
69 StoreObjPtr<EventMetaData> eventMetaData("EventMetaData", DataStore::c_Event);
70 Float_t event_num = eventMetaData->getEvent();
71 Float_t event_run = eventMetaData->getRun();
72 Float_t event_exp = eventMetaData->getExperiment();
73 Float_t event_prod = eventMetaData->getProduction();
74
75 B2DEBUG(29, "+++++ Loop on Tracks");
77
78 for (Track& track : tracks) {
79
80 RecoTrack* recoTrack = track.getRelationsTo<RecoTrack>()[0];
81 if (recoTrack == nullptr) {
82 // if no recoTrack is associated to Track, we skip the track
83 B2WARNING(" the RecoTrack associated to Track is nullptr!");
84 continue;
85 }
86
87 Float_t nhits = recoTrack->getNumberOfTrackingHits();
88 Float_t ncdc = recoTrack->getNumberOfCDCHits();
89 Float_t npxd = recoTrack->getNumberOfPXDHits();
90 Float_t nsvd = recoTrack->getNumberOfSVDHits();
91
92 const Const::ChargedStable pdg_list[4] = {Const::pion, Const::kaon, Const::proton, Const::deuteron}; // loop only on these hypotheses
93
94 std::map <Const::ChargedStable, Float_t> flag; // is the particle hypothesis existing in the track?
95 std::map <Const::ChargedStable, Float_t> trk_x, trk_y, trk_z, trk_px, trk_py, trk_pz, trk_p, trk_pt, trk_theta, trk_phi;
96 std::map <Const::ChargedStable, Float_t> trk_charge, trk_chi2, trk_ndf, trk_pvalue, trk_nfailed;
97 std::map <Const::ChargedStable, Float_t> nhits_pid, ncdc_pid, nsvd_pid, npxd_pid;
98 std::map <Const::ChargedStable, Float_t> first_cdc, last_cdc, first_svd, last_svd;
99
100 for (const Const::ChargedStable& pdgIter : pdg_list) {
101 trk_x[pdgIter] = 0.; trk_y[pdgIter] = 0.; trk_z[pdgIter] = 0.;
102 trk_px[pdgIter] = 0.; trk_py[pdgIter] = 0.; trk_pz[pdgIter] = 0.;
103 trk_p[pdgIter] = 0.; trk_pt[pdgIter] = 0.; trk_theta[pdgIter] = 0.; trk_phi[pdgIter] = 0.;
104 trk_charge[pdgIter] = 0.;
105 trk_chi2[pdgIter] = 0.; trk_ndf[pdgIter] = 0.; trk_pvalue[pdgIter] = 0.; trk_nfailed[pdgIter] = 0.;
106 nhits_pid[pdgIter] = 0.; ncdc_pid[pdgIter] = 0.; nsvd_pid[pdgIter] = 0.; npxd_pid[pdgIter] = 0.;
107 first_cdc[pdgIter] = -100.; last_cdc [pdgIter] = -100.;
108 first_svd[pdgIter] = -100.; last_svd [pdgIter] = -100.;
109
110 const TrackFitResult* fitResult = track.getTrackFitResult(pdgIter);
111 if ((fitResult != nullptr) && (fitResult->getParticleType() == pdgIter)) {
112 flag[pdgIter] = kTRUE;
113 } else {
114 flag[pdgIter] = kFALSE;
115 continue;
116 }
117
118 trk_x[pdgIter] = fitResult->getPosition().X();
119 trk_y[pdgIter] = fitResult->getPosition().Y();
120 trk_z[pdgIter] = fitResult->getPosition().Z();
121 trk_px[pdgIter] = fitResult->getMomentum().X();
122 trk_py[pdgIter] = fitResult->getMomentum().Y();
123 trk_pz[pdgIter] = fitResult->getMomentum().Z();
124 trk_p[pdgIter] = fitResult->getMomentum().R();
125 trk_pt[pdgIter] = fitResult->getMomentum().Rho();
126 trk_theta[pdgIter] = fitResult->getMomentum().Theta() * TMath::RadToDeg();
127 trk_phi[pdgIter] = fitResult->getMomentum().Phi() * TMath::RadToDeg();
128 trk_charge[pdgIter] = fitResult->getChargeSign();
129 double chi2 = recoTrack->getTrackFitStatus(recoTrack->getTrackRepresentationForPDG(pdgIter.getPDGCode()))->getChi2();
130 if (isnan(chi2)) chi2 = -10;
131 if (isinf(chi2)) chi2 = -20;
132 trk_chi2[pdgIter] = chi2;
133 trk_ndf[pdgIter] = fitResult->getNDF();
134 trk_pvalue[pdgIter] = fitResult->getPValue();
135 trk_nfailed[pdgIter] = recoTrack->getTrackFitStatus(recoTrack->getTrackRepresentationForPDG(
136 pdgIter.getPDGCode()))->getNFailedPoints();
137 ncdc_pid[pdgIter] = fitResult->getHitPatternCDC().getNHits();
138 npxd_pid[pdgIter] = fitResult->getHitPatternVXD().getNPXDHits();
139 nsvd_pid[pdgIter] = fitResult->getHitPatternVXD().getNSVDHits();
140 nhits_pid[pdgIter] = ncdc_pid[pdgIter] + npxd_pid[pdgIter] + nsvd_pid[pdgIter];
141 first_cdc[pdgIter] = fitResult->getHitPatternCDC().getFirstLayer();
142 last_cdc[pdgIter] = fitResult->getHitPatternCDC().getLastLayer();
143 first_svd[pdgIter] = fitResult->getHitPatternVXD().getFirstSVDLayer();
144 last_svd[pdgIter] = fitResult->getHitPatternVXD().getLastSVDLayer();
145 }
146
147 Float_t buffer[] = {event_num, event_run, event_exp, event_prod,
148 nhits, ncdc, npxd, nsvd,
149 (Float_t)recoTrack->getPositionSeed().X(), (Float_t)recoTrack->getPositionSeed().Y(), (Float_t)recoTrack->getPositionSeed().Z(),
150 (Float_t)recoTrack->getMomentumSeed().X(), (Float_t)recoTrack->getMomentumSeed().Y(), (Float_t)recoTrack->getMomentumSeed().Z(), (Float_t)recoTrack->getMomentumSeed().R(), (Float_t)recoTrack->getMomentumSeed().Rho(),
151 (Float_t)(recoTrack->getMomentumSeed().Theta()* TMath::RadToDeg()), (Float_t)(recoTrack->getMomentumSeed().Phi()* TMath::RadToDeg()), (Float_t)recoTrack->getChargeSeed(),
152 nhits_pid[Const::pion], ncdc_pid[Const::pion], npxd_pid[Const::pion], nsvd_pid[Const::pion],
153 nhits_pid[Const::kaon], ncdc_pid[Const::kaon], npxd_pid[Const::kaon], nsvd_pid[Const::kaon],
154 nhits_pid[Const::proton], ncdc_pid[Const::proton], npxd_pid[Const::proton], nsvd_pid[Const::proton],
155 nhits_pid[Const::deuteron], ncdc_pid[Const::deuteron], npxd_pid[Const::deuteron], nsvd_pid[Const::deuteron],
156 flag[Const::pion], flag[Const::kaon], flag[Const::proton], flag[Const::deuteron],
157 trk_x[Const::pion], trk_y[Const::pion], trk_z[Const::pion],
158 trk_px[Const::pion], trk_py[Const::pion], trk_pz[Const::pion], trk_p[Const::pion], trk_pt[Const::pion], trk_theta[Const::pion],
159 trk_phi[Const::pion], trk_charge[Const::pion], trk_chi2[Const::pion], trk_ndf[Const::pion], trk_pvalue[Const::pion], trk_nfailed[Const::pion],
160 trk_x[Const::kaon], trk_y[Const::kaon], trk_z[Const::kaon],
161 trk_px[Const::kaon], trk_py[Const::kaon], trk_pz[Const::kaon], trk_p[Const::kaon], trk_pt[Const::kaon], trk_theta[Const::kaon],
162 trk_phi[Const::kaon], trk_charge[Const::kaon], trk_chi2[Const::kaon], trk_ndf[Const::kaon], trk_pvalue[Const::kaon], trk_nfailed[Const::kaon],
163 trk_x[Const::proton], trk_y[Const::proton], trk_z[Const::proton],
164 trk_px[Const::proton], trk_py[Const::proton], trk_pz[Const::proton], trk_p[Const::proton], trk_pt[Const::proton], trk_theta[Const::proton],
165 trk_phi[Const::proton], trk_charge[Const::proton], trk_chi2[Const::proton], trk_ndf[Const::proton], trk_pvalue[Const::proton], trk_nfailed[Const::proton],
166 trk_x[Const::deuteron], trk_y[Const::deuteron], trk_z[Const::deuteron],
167 trk_px[Const::deuteron], trk_py[Const::deuteron], trk_pz[Const::deuteron], trk_p[Const::deuteron], trk_pt[Const::deuteron], trk_theta[Const::deuteron],
168 trk_phi[Const::deuteron], trk_charge[Const::deuteron], trk_chi2[Const::deuteron], trk_ndf[Const::deuteron], trk_pvalue[Const::deuteron], trk_nfailed[Const::deuteron],
169 first_cdc[Const::pion], last_cdc[Const::pion], first_svd[Const::pion], last_svd[Const::pion],
170 first_cdc[Const::kaon], last_cdc[Const::kaon], first_svd[Const::kaon], last_svd[Const::kaon],
171 first_cdc[Const::proton], last_cdc[Const::proton], first_svd[Const::proton], last_svd[Const::proton],
172 first_cdc[Const::deuteron], last_cdc[Const::deuteron], first_svd[Const::deuteron], last_svd[Const::deuteron]
173 };
174 m_n_MultiParticle->Fill(buffer);
175 }
176
177
178
179}
180
181
183{
184
185 if (m_rootFilePtr != nullptr) {
186 m_rootFilePtr->cd();
187 m_n_MultiParticle->Write();
188
189 m_rootFilePtr->Close();
190 }
191
192}
Provides a type-safe way to pass members of the chargedStableSet set.
Definition: Const.h:589
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 deuteron
deuteron particle
Definition: Const.h:664
@ c_Event
Different object in each event, all objects/arrays are invalidated after event() function has been ca...
Definition: DataStore.h:59
int m_ParticleHypothesis
Particle Hypothesis for the track fit (default: 211)
std::string m_TracksName
Track StoreArray name.
void initialize() override
Require the store arrays and create the output root file.
void event() override
Loop over Track objects and fill ntuples with tracking parameters.
void terminate() override
Save output root file with ntuple.
TNtuple * m_n_MultiParticle
Multi particle ntuple.
std::string m_RecoTracksName
RecoTrack StoreArray name.
FillTrackFitNtupleModule()
Default Empty Constructor.
StoreArray< RecoTrack > m_RecoTracks
RecoTrack StoreArray.
StoreArray< Track > m_Tracks
Track StoreArray.
short getLastLayer() const
Returns the index of the last layer with a hit.
unsigned short getNHits() const
Get the total Number of CDC hits in the fit.
short getFirstLayer() const
Returns the index of the first layer with a hit.
short getLastSVDLayer() const
Get the last activated SVD layer index.
short getFirstSVDLayer() const
Get the first activated SVD layer index.
unsigned short getNPXDHits() const
Get total number of hits in the PXD.
unsigned short getNSVDHits() const
Get total number of hits in the SVD.
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
TFile * m_rootFilePtr
pointer at root file used for storing histograms
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
ROOT::Math::XYZVector getPositionSeed() const
Return the position seed stored in the reco track. ATTENTION: This is not the fitted position.
Definition: RecoTrack.h:480
unsigned int getNumberOfSVDHits() const
Return the number of svd hits.
Definition: RecoTrack.h:424
unsigned int getNumberOfCDCHits() const
Return the number of cdc hits.
Definition: RecoTrack.h:427
unsigned int getNumberOfTrackingHits() const
Return the number of cdc + svd + pxd hits.
Definition: RecoTrack.h:443
short int getChargeSeed() const
Return the charge seed stored in the reco track. ATTENTION: This is not the fitted charge.
Definition: RecoTrack.h:508
ROOT::Math::XYZVector getMomentumSeed() const
Return the momentum seed stored in the reco track. ATTENTION: This is not the fitted momentum.
Definition: RecoTrack.h:487
genfit::AbsTrackRep * getTrackRepresentationForPDG(int pdgCode) const
Return an already created track representation of the given reco track for the PDG.
Definition: RecoTrack.cc:475
unsigned int getNumberOfPXDHits() const
Return the number of pxd hits.
Definition: RecoTrack.h:421
const genfit::FitStatus * getTrackFitStatus(const genfit::AbsTrackRep *representation=nullptr) const
Return the track fit status for the given representation or for the cardinal one. You are not allowed...
Definition: RecoTrack.h:621
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
Values of the result of a track fit with a given particle hypothesis.
float getNDF() const
Getter for number of degrees of freedom of the track fit.
short getChargeSign() const
Return track charge (1 or -1).
double getPValue() const
Getter for Chi2 Probability of the track fit.
Const::ParticleType getParticleType() const
Getter for ParticleType of the mass hypothesis of the track fit.
ROOT::Math::XYZVector getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
ROOT::Math::XYZVector getPosition() const
Getter for vector of position at closest approach of track in r/phi projection.
HitPatternCDC getHitPatternCDC() const
Getter for the hit pattern in the CDC;.
HitPatternVXD getHitPatternVXD() const
Getter for the hit pattern in the VXD;.
Class that bundles various TrackFitResults.
Definition: Track.h:25
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
Abstract base class for different kinds of events.