Belle II Software  release-06-00-14
V0findingPerformanceEvaluationModule.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/V0findingPerformanceEvaluationModule.h>
10 
11 #include <tracking/dataobjects/MCParticleInfo.h>
12 #include <tracking/dataobjects/V0ValidationVertex.h>
13 
14 #include <framework/datastore/StoreArray.h>
15 #include <framework/datastore/RelationVector.h>
16 
17 #include <mdst/dataobjects/Track.h>
18 
19 #include <framework/gearbox/Const.h>
20 
21 #include <framework/geometry/BFieldManager.h>
22 
23 #include <root/TAxis.h>
24 #include <root/TObject.h>
25 
26 #include <boost/foreach.hpp>
27 #include <vector>
28 
29 using namespace Belle2;
30 
31 //-----------------------------------------------------------------
32 // Register the Module
33 //-----------------------------------------------------------------
34 REG_MODULE(V0findingPerformanceEvaluation)
35 
36 V0findingPerformanceEvaluationModule::V0findingPerformanceEvaluationModule() :
37  Module()
38 {
39 
40  setDescription("This module evaluates the V0 finding package performance");
41 
42  addParam("outputFileName", m_rootFileName, "Name of output root file.",
43  std::string("V0findingPerformanceEvaluation_output.root"));
44  addParam("V0sName", m_V0sName, "Name of V0 collection.", std::string("V0ValidationVertexs"));
45  addParam("MCParticlesName", m_MCParticlesName, "Name of MC Particle collection.", std::string(""));
46 
47 }
48 
49 V0findingPerformanceEvaluationModule::~V0findingPerformanceEvaluationModule()
50 {
51 
52 }
53 
55 {
56  StoreArray<MCParticle> mcParticles;
57  mcParticles.isRequired(m_MCParticlesName);
58 
59  StoreArray<V0ValidationVertex> v0ValidationVertices;
60  v0ValidationVertices.isRequired(m_V0sName);
61 
62  //create list of histograms to be saved in the rootfile
63  m_histoList = new TList;
64  m_histoList_multiplicity = new TList;
65  m_histoList_efficiency = new TList;
66  m_histoList_purity = new TList;
67  m_histoList_trkQuality = new TList;
68 
69  //set the ROOT File
70  m_rootFilePtr = new TFile(m_rootFileName.c_str(), "RECREATE");
71 
72  //now create histograms
73 
74  //multiplicity histograms
75  m_multiplicityV0s = createHistogram1D("h1nV0", "number of V0s per MC Particle", 8, -0.5, 7.5, "# V0s", m_histoList_multiplicity);
76 
77  m_multiplicityMCParticles = createHistogram1D("h1nMCPrtcl", "number of MCParticles per V0s", 5, -0.5, 4.5,
78  "# MCParticles", m_histoList_multiplicity);
79 
80  m_MCParticlePDGcode = createHistogram1D("h1PDGcode", "PDG code of MCParticles", 6244, Const::antiLambda.getPDGCode(),
81  Const::Lambda.getPDGCode(),
82  "PDG code", m_histoList_multiplicity);
83 
84 
85  //vertex and momentum parameters errors
86  m_h1_vtxX_err = createHistogram1D("h1vtxXerr", "vtxX error", 100, 0, 0.1, "#sigma_{vtxX} (cm)", m_histoList);
87  m_h1_vtxY_err = createHistogram1D("h1vtxYerr", "vtxY error", 100, 0, 0.1, "#sigma_{vtxY} (cm)", m_histoList);
88  m_h1_vtxZ_err = createHistogram1D("h1vtxZerr", "vtxZ error", 100, 0, 0.3, "#sigma_{vtxZ} (cm)", m_histoList);
89  m_h2_vtxTvsR_err = createHistogram2D("h2vtxTerrVsR", "vtxT error vs R", 100, 0, 100, "R (cm)", 100, 0, 0.3, "#sigma_{vtxT} (cm)",
90  m_histoList);
91  // m_h1_mom_err = createHistogram1D("h1momerr", "mom error", 100, 0, 0.1, "#sigma_{p} (GeV/c)", m_histoList);
92  // m_h1_mass_err = createHistogram1D("h1masserr", "mass error", 100, 0, 1, "#sigma_{m} (GeV/c2)", m_histoList);
93  //vertex and momentum parameters residuals
94  m_h1_vtxX_res = createHistogram1D("h1vtxXres", "vtxX resid", 100, -0.2, 0.2, "vtxX resid (cm)", m_histoList);
95  m_h1_vtxY_res = createHistogram1D("h1vtxYres", "vtxY resid", 100, -0.2, 0.2, "vtxY resid (cm)", m_histoList);
96  m_h1_vtxZ_res = createHistogram1D("h1vtxZres", "vtxZ resid", 100, -0.5, 0.5, "vtxZ resid (cm)", m_histoList);
97  m_h1_mom_res = createHistogram1D("h1momres", "mom resid", 1000, -0.5, 0.5, "mom resid (GeV/c)", m_histoList);
98  m_h1_mass_res = createHistogram1D("h1massres", "mass resid", 500, -0.3, 0.3, "mass resid (GeV/c)", m_histoList);
99 
100  //vertex and momentum parameters pulls
101  m_h1_vtxX_pll = createHistogram1D("h1vtxXpll", "vtxX pull", 100, -5, 5, "vtxX pull", m_histoList);
102  m_h1_vtxY_pll = createHistogram1D("h1vtxYpll", "vtxY pull", 100, -5, 5, "vtxY pull", m_histoList);
103  m_h1_vtxZ_pll = createHistogram1D("h1vtxZpll", "vtxZ pull", 100, -5, 5, "vtxZ pull", m_histoList);
104  // m_h1_mom_pll = createHistogram1D("h1mompll", "mom pull", 100, -5, 5, "momX pull", m_histoList);
105  // m_h1_mass_pll = createHistogram1D("h1masspll", "mass pull", 100, -5, 5, "momY pull", m_histoList);
106 
107 
108  m_h1_ChiSquare = createHistogram1D("h1Chi2", "Chi2 of the fit", 100, 0, 20, "Chi2", m_histoList_trkQuality);
109 
110  m_h1_nMatchedDau = createHistogram1D("h1nMatchedDau", "Number of Matched MCParticle Daughters", 3, -0.5, 2.5, "# matched dau",
111  m_histoList);
112 
113 
114  m_h2_mom = createHistogram2D("h2mom", "reco VS true momentum", 100, 0, 3, "V0 mom (GeV/c)", 100, 0, 3, "MC mom (GeV/c)",
115  m_histoList);
116  m_h2_mass = createHistogram2D("h2mass", "reco VS true mass", 100, 0, 1.5, "V0 mass (GeV/c2)", 100, 0, 1.5, "MC mass (GeV/c)",
117  m_histoList);
118 
119  //histograms to produce efficiency plots
120  Double_t bins_pt[9 + 1] = {0, 0.05, 0.1, 0.15, 0.2, 0.3, 0.5, 1, 2, 3.5}; //GeV/c
121  Double_t bins_theta[10 + 1] = {0, 0.25, 0.5, 0.75, 0.75 + 0.32, 0.75 + 2 * 0.32, 0.75 + 3 * 0.32, 0.75 + 4 * 0.32, 0.75 + 5 * 0.32, 2.65, TMath::Pi()};
122  Double_t bins_phi[14 + 1];
123  Double_t width_phi = 2 * TMath::Pi() / 14;
124  for (int bin = 0; bin < 14 + 1; bin++)
125  bins_phi[bin] = - TMath::Pi() + bin * width_phi;
126 
127  m_h1_MCParticle_R = createHistogram1D("h1nMCParticleVSr", "entry per MCParticles", 50, 0, 20, "transverse L", m_histoList);
128 
129  m_h1_V0sPerMCParticle_R = (TH1F*)duplicateHistogram("h1nV0perMCvsR", "entry per V0 related to a MCParticle", m_h1_MCParticle_R,
130  m_histoList);
131 
132 
133  m_h3_MCParticle = createHistogram3D("h3MCParticle", "entry per MCParticle",
134  9, bins_pt, "p_{t} (GeV/c)",
135  10, bins_theta, "#theta",
136  14, bins_phi, "#phi" /*, m_histoList*/);
137 
138  m_h3_V0sPerMCParticle = (TH3F*)duplicateHistogram("h3V0sPerMCParticle",
139  "entry per V0 connected to a MCParticle",
140  m_h3_MCParticle /*, m_histoList*/);
141 
142  m_h3_V0s = (TH3F*)duplicateHistogram("h3V0s", "entry per V0",
143  m_h3_MCParticle /*, m_histoList*/);
144 
145  //histograms to produce purity plots
146  m_h3_MCParticlesPerV0 = (TH3F*)duplicateHistogram("h3MCParticlesPerV0",
147  "entry per MCParticle connected to a V0",
148  m_h3_MCParticle /*, m_histoList*/);
149 }
150 
152 {
153 
154 }
155 
157 {
158 
160 
161  B2Vector3D magField = BFieldManager::getField(0, 0, 0) / Unit::T;
162 
163  B2DEBUG(99, "+++++ 1. loop on MCParticles");
164  BOOST_FOREACH(MCParticle & mcParticle, mcParticles) {
165 
166  if (! isV0(mcParticle))
167  continue;
168 
169  int nMatchedDau = nMatchedDaughters(mcParticle);
170  m_h1_nMatchedDau->Fill(nMatchedDau);
171 
172  //proceed only in case the MCParticle daughters have both one associated reconstructed track:
173  if (nMatchedDau != 2)
174  continue;
175 
176  int pdgCode = mcParticle.getPDG();
177  B2DEBUG(99, "MCParticle has PDG code " << pdgCode);
178  m_MCParticlePDGcode->Fill(mcParticle.getPDG());
179 
180  MCParticleInfo mcParticleInfo(mcParticle, magField);
181 
182  TVector3 MC_prodvtx = mcParticle.getVertex();
183  TVector3 MC_vtx = mcParticle.getDecayVertex();
184  float MC_mom = mcParticle.getMomentum().Mag();
185  float MC_mass = mcParticle.getMass();
186  TVector3 MC_FL = MC_vtx - MC_prodvtx;
187  float flightR = sqrt(MC_FL.X() * MC_FL.X() + MC_FL.Y() * MC_FL.Y());
188  m_h1_MCParticle_R->Fill(flightR);
189 
190  m_h3_MCParticle->Fill(mcParticleInfo.getPt(), mcParticleInfo.getPtheta(), mcParticleInfo.getPphi());
191 
192  //1. retrieve all the V0s related to the MCParticle
193 
194  //1.0 check if there is a V0
195  RelationVector<V0ValidationVertex> V0s_toMCParticle =
196  DataStore::getRelationsWithObj<V0ValidationVertex>(&mcParticle, m_V0sName);
197 
198  m_multiplicityV0s->Fill(V0s_toMCParticle.size());
199 
200  if (V0s_toMCParticle.size() > 0)
201  m_h3_V0sPerMCParticle->Fill(mcParticleInfo.getPt(), mcParticleInfo.getPtheta(), mcParticleInfo.getPphi());
202 
203  for (int v0 = 0; v0 < (int)V0s_toMCParticle.size(); v0++) {
204 
205  TVector3 V0_vtx = V0s_toMCParticle[v0]->getVertexPosition();
206  float V0_mom = V0s_toMCParticle[v0]->getFittedMomentum();
207  float V0_chi2 = V0s_toMCParticle[v0]->getVertexChi2();
208  float V0_mass = V0s_toMCParticle[v0]->getFittedInvariantMass();
209  TMatrixDSym V0_cov = V0s_toMCParticle[v0]->getVertexPositionCovariance();
210 
211  m_h1_vtxX_err->Fill(sqrt(V0_cov[0][0]));
212  m_h1_vtxY_err->Fill(sqrt(V0_cov[1][1]));
213  m_h1_vtxZ_err->Fill(sqrt(V0_cov[2][2]));
214  m_h2_vtxTvsR_err->Fill(flightR, sqrt(V0_cov[0][0] + V0_cov[1][1]));
215 
216  m_h1_V0sPerMCParticle_R->Fill(flightR);
217 
218  m_h1_vtxX_res->Fill(V0_vtx.X() - MC_vtx.X());
219  m_h1_vtxY_res->Fill(V0_vtx.Y() - MC_vtx.Y());
220  m_h1_vtxZ_res->Fill(V0_vtx.Z() - MC_vtx.Z());
221 
222  m_h1_mom_res->Fill(V0_mom - MC_mom);
223  m_h2_mom->Fill(V0_mom, MC_mom);
224  m_h1_mass_res->Fill(V0_mass - MC_mass);
225  m_h2_mass->Fill(V0_mass, MC_mass);
226 
227  m_h1_vtxX_pll->Fill((V0_vtx.X() - MC_vtx.X()) / sqrt(V0_cov[0][0]));
228  m_h1_vtxY_pll->Fill((V0_vtx.Y() - MC_vtx.Y()) / sqrt(V0_cov[1][1]));
229  m_h1_vtxZ_pll->Fill((V0_vtx.Z() - MC_vtx.Z()) / sqrt(V0_cov[2][2]));
230 
231  m_h1_ChiSquare->Fill(V0_chi2);
232 
233  }
234 
235 
236  }
237 
238 
239  B2DEBUG(99, "+++++ 2. loop on V0s");
240 
241 
243 
244 
245  BOOST_FOREACH(V0ValidationVertex & v0, V0s) {
246 
247  int nMCParticles = 0;
248 
249  // m_h3_V0s>Fill(mcParticleInfo.getPt(), mcParticleInfo.getPtheta(), mcParticleInfo.getPphi());
250  //2. retrieve all the MCParticles related to the V0s
251  RelationVector<MCParticle> MCParticles_fromV0 =
252  DataStore::getRelationsWithObj<MCParticle>(&v0, m_MCParticlesName);
253 
254  nMCParticles = MCParticles_fromV0.size();
255 
256  if (nMCParticles == 0)
257  continue;
258 
259  MCParticleInfo mcParticleInfo(* MCParticles_fromV0[0], magField);
260  m_h3_MCParticlesPerV0->Fill(mcParticleInfo.getPt(), mcParticleInfo.getPtheta(), mcParticleInfo.getPphi());
261  m_multiplicityMCParticles->Fill(nMCParticles);
262 
263  }
264 
265 }
266 
268 {
269 
270  double num = 0;
271  double den = 0;
272 
273  for (int bin = 1; bin < m_multiplicityV0s->GetNbinsX(); bin ++)
274  num += m_multiplicityV0s->GetBinContent(bin + 1);
275  den = m_multiplicityV0s->GetEntries();
276  double efficiency = num / den ;
277  double efficiencyErr = sqrt(efficiency * (1 - efficiency)) / sqrt(den);
278 
279  double nMCParticles = 0;
280  for (int bin = 1; bin < m_multiplicityMCParticles->GetNbinsX(); bin ++)
281  nMCParticles += m_multiplicityMCParticles->GetBinContent(bin + 1);
282  double purity = nMCParticles / m_multiplicityMCParticles->GetEntries();
283  double purityErr = sqrt(purity * (1 - purity)) / sqrt(m_multiplicityMCParticles->GetEntries());
284 
285  B2INFO("");
286  B2INFO("~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~");
287  B2INFO("~ V0 Finding Performance Evaluation ~ SHORT SUMMARY ~");
288  B2INFO("");
289  B2INFO(" + overall:");
290  B2INFO(" efficiency = (" << efficiency * 100 << " +/- " << efficiencyErr * 100 << ")% ");
291  B2INFO(" purity = (" << purity * 100 << " +/- " << purityErr * 100 << ")% ");
292  B2INFO("");
293  B2INFO("~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~");
294 }
295 
297 {
298 
299  TH1F* h_eff_R = (TH1F*)duplicateHistogram("h_eff_R", "efficiency vs R", m_h1_MCParticle_R, m_histoList_efficiency);
300 
301  for (int bin = 0; bin < h_eff_R->GetXaxis()->GetNbins(); bin++) {
302  float num = m_h1_V0sPerMCParticle_R->GetBinContent(bin + 1);
303  float den = m_h1_MCParticle_R->GetBinContent(bin + 1);
304  double eff = 0;
305  double err = 0;
306 
307  if (den > 0) {
308  eff = (double)num / den;
309  err = sqrt(eff * (1 - eff)) / sqrt(den);
310  }
311 
312  h_eff_R->SetBinContent(bin + 1, eff);
313  h_eff_R->SetBinError(bin + 1, err);
314  }
315 
317 
319 
320  // addPurityPlots(m_histoList_purity, m_h3_MCParticlesPerV0, m_h3_V0s);
321 
322  if (m_rootFilePtr != nullptr) {
323  m_rootFilePtr->cd();
324 
325  TDirectory* oldDir = gDirectory;
326 
327  TDirectory* dir_multiplicity = oldDir->mkdir("multiplicity");
328  dir_multiplicity->cd();
329  TIter nextH_multiplicity(m_histoList_multiplicity);
330  TObject* obj;
331  while ((obj = nextH_multiplicity()))
332  obj->Write();
333 
334  TDirectory* dir_efficiency = oldDir->mkdir("efficiency");
335  dir_efficiency->cd();
336  TIter nextH_efficiency(m_histoList_efficiency);
337  while ((obj = nextH_efficiency()))
338  obj->Write();
339 
340  TDirectory* dir_purity = oldDir->mkdir("purity");
341  dir_purity->cd();
342  TIter nextH_purity(m_histoList_purity);
343  while ((obj = nextH_purity()))
344  obj->Write();
345 
346  TDirectory* dir_trkQuality = oldDir->mkdir("trkQuality");
347  dir_trkQuality->cd();
348  TIter nextH_trkQuality(m_histoList_trkQuality);
349  while ((obj = nextH_trkQuality()))
350  obj->Write();
351 
352  TIter nextH(m_histoList);
353  while ((obj = nextH()))
354  obj->Write();
355 
356 
357  m_rootFilePtr->Close();
358  }
359 
360 }
361 
362 
364 {
365 
366  bool isGamma = false;
367  if (abs(the_mcParticle.getPDG()) == Const::photon.getPDGCode())
368  isGamma = true;
369 
370  bool isK_S0 = false;
371  if (abs(the_mcParticle.getPDG()) == Const::Kshort.getPDGCode())
372  isK_S0 = true;
373 
374  bool isK_0 = false;
375  if (abs(the_mcParticle.getPDG()) == 311)
376  isK_0 = true;
377 
378  bool isLambda = false;
379  if (abs(the_mcParticle.getPDG()) == Const::Lambda.getPDGCode())
380  isLambda = true;
381 
382  bool twoProngs = false;
383  bool twoChargedProngs = false;
384 
385  if (the_mcParticle.getDaughters().size() == 2)
386  twoProngs = true;
387 
388  if (twoProngs)
389  if (the_mcParticle.getDaughters()[0]->getCharge() * the_mcParticle.getDaughters()[1]->getCharge() < 0)
390  twoChargedProngs = true;
391 
392  return ((isGamma || isK_S0 || isK_0 || isLambda) && twoChargedProngs);
393 
394 }
395 
397 {
398 
399  int nMatchedDau = 0;
400 
401  std::vector< MCParticle* > MCPart_dau = the_mcParticle.getDaughters();
402 
403  bool first = false;
404  bool second = false;
405 
406  RelationVector<Track> Tracks_fromMCParticle_0 = DataStore::getRelationsWithObj<Track>(MCPart_dau[0]);
407  if (Tracks_fromMCParticle_0.size() > 0)
408  first = true;
409 
410  RelationVector<Track> Tracks_fromMCParticle_1 = DataStore::getRelationsWithObj<Track>(MCPart_dau[1]);
411  if (Tracks_fromMCParticle_1.size() > 0)
412  second = true;
413 
414 
415  if (first)
416  nMatchedDau++;
417 
418  if (second)
419  nMatchedDau++;
420 
421 
422  return nMatchedDau;
423 
424 }
int getPDGCode() const
PDG code.
Definition: Const.h:354
static const ParticleType Lambda
Lambda particle.
Definition: Const.h:559
static const ParticleType antiLambda
Anti-Lambda particle.
Definition: Const.h:560
static const ParticleType Kshort
K^0_S particle.
Definition: Const.h:557
static const ParticleType photon
photon particle
Definition: Const.h:554
This struct is used by the TrackingPerformanceEvaluation Module to save information of reconstructed ...
double getPt()
Getter for transverse momentum.
double getPtheta()
Getter for theta of momentum vector.
double getPphi()
Getter for phi of momentum vector.
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
std::vector< Belle2::MCParticle * > getDaughters() const
Get vector of all daughter particles, empty vector if none.
Definition: MCParticle.cc:50
float getMass() const
Return the particle mass in GeV.
Definition: MCParticle.h:135
TVector3 getVertex() const
Return production vertex position, shorthand for getProductionVertex().
Definition: MCParticle.h:183
TVector3 getMomentum() const
Return momentum.
Definition: MCParticle.h:198
TVector3 getDecayVertex() const
Return decay vertex.
Definition: MCParticle.h:219
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
Base class for Modules.
Definition: Module.h:72
TH1F * createHistogram1D(const char *name, const char *title, Int_t nbins, Double_t min, Double_t max, const char *xtitle, TList *histoList=nullptr)
Create a 1D histogram and add it to the TList of 1D-histograms.
TH1 * duplicateHistogram(const char *newname, const char *newtitle, TH1 *h, TList *histoList=nullptr)
Make a copy of a 1D histogram and add it to the TList of 1D-histograms.
void addInefficiencyPlots(TList *graphList=nullptr, TH3F *h3_xPerMCParticle=nullptr, TH3F *h3_MCParticle=nullptr)
Create pt-, theta- and phi-inefficiency 1D histograms and add them to the TList of 1D-histograms.
TList * m_histoList_trkQuality
List of track-quality histograms.
TList * m_histoList_purity
List of purity histograms.
TList * m_histoList
List of performance-evaluation histograms.
TList * m_histoList_multiplicity
List of multiplicity histograms.
void addEfficiencyPlots(TList *graphList=nullptr, TH3F *h3_xPerMCParticle=nullptr, TH3F *h3_MCParticle=nullptr)
Create pt-, theta- and phi-efficiency 1D histograms and add them to the TList of 1D-histograms.
TH3F * createHistogram3D(const char *name, const char *title, Int_t nbinsX, Double_t minX, Double_t maxX, const char *titleX, Int_t nbinsY, Double_t minY, Double_t maxY, const char *titleY, Int_t nbinsZ, Double_t minZ, Double_t maxZ, const char *titleZ, TList *histoList=nullptr)
Create a 3D histogram and add it to the TList of 3D-histograms.
TH2F * createHistogram2D(const char *name, const char *title, Int_t nbinsX, Double_t minX, Double_t maxX, const char *titleX, Int_t nbinsY, Double_t minY, Double_t maxY, const char *titleY, TList *histoList=nullptr)
Create a 2D histogram and add it to the TList of 2D-histograms.
TFile * m_rootFilePtr
pointer at root file used for storing histograms
TList * m_histoList_efficiency
List of efficiency histograms.
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
static const double T
[tesla]
Definition: Unit.h:120
Class which stores some additional information on V0 vertices.
TH2F * m_h2_vtxTvsR_err
vtx error on transverse plane VS transverse flight length
void event() override
This method is the core of the module.
void endRun() override
This method is called if the current run ends.
void terminate() override
This method is called at the end of the event processing.
int nMatchedDaughters(const MCParticle &the_mcParticle)
number of truth matched dauhters
void beginRun() override
Called when entering a new run.
TH1F * m_multiplicityMCParticles
number of MCParticles per fitted Track
#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.
Abstract base class for different kinds of events.