Belle II Software  release-05-02-19
MCTrackCandClassifierModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2011 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Giulia Casarosa *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <tracking/modules/mcTrackCandClassifier/MCTrackCandClassifierModule.h>
12 
13 #include <pxd/dataobjects/PXDTrueHit.h>
14 #include <pxd/dataobjects/PXDCluster.h>
15 #include <svd/dataobjects/SVDTrueHit.h>
16 #include <svd/dataobjects/SVDCluster.h>
17 #include <vxd/geometry/GeoCache.h>
18 
19 #include <framework/datastore/StoreArray.h>
20 #include <framework/geometry/B2Vector3.h>
21 #include <framework/geometry/BFieldManager.h>
22 
23 #include <genfit/TrackCand.h>
24 
25 #include <boost/foreach.hpp>
26 
27 #include <TH2F.h>
28 
29 using namespace std;
30 using namespace Belle2;
31 
33 REG_MODULE(MCTrackCandClassifier)
34 
35 MCTrackCandClassifierModule::MCTrackCandClassifierModule() : Module()
37  , m_rootFilePtr(NULL)
38 {
39  //Set module properties
40  setDescription("This module is meant to classify the MCTrackCands as either ideal, fine and nasty");
42 
43  addParam("MCParticlesName", m_mcParticlesName,
44  "Name of MC Particle collection.",
45  std::string(""));
46 
47  addParam("MCTrackCandCollName", m_mcTrackCandsColName,
48  "Name of the input collection of MC track candidates",
49  std::string(""));
50 
51  addParam("rootFileName", m_rootFileName,
52  "Name of the root file",
53  std::string("MCTrackCandClassifier.root"));
54 
55  addParam("isInAnnulusCriterium", m_applyAnnulus,
56  "Require that the hit is in the expected annulus",
57  bool(true));
58  addParam("isInSemiplaneCriterium", m_applySemiplane,
59  "Require that the hit is in the expected semiplane",
60  bool(true));
61  addParam("isInFirstLapCriterium", m_applyLap,
62  "Require that the hit belong to the first lap in the transverse plane",
63  bool(true));
64  addParam("isInWedgePartCriterium", m_applyWedge,
65  "Require that the hit belong to the barrel part of the SVD",
66  bool(true));
67  addParam("removeBadHits", m_removeBadHits,
68  "Remove the clusters that do not satisfy the criteria from the idealMCTrackCands",
69  bool(true));
70 
71  addParam("minNhits", m_minHit,
72  "Minimum number of 1D Clusters to classify the MCTrackCand as ideal",
73  int(5));
74 
75  addParam("nSigma_dR", m_nSigma, "n sigma dR", int(3));
76 
77  addParam("lapFraction", m_fraction, "Fraction of lap", double(1));
78 
79  addParam("usePXD", m_usePXD, "Use the PXD or not", bool(true));
80 }
81 
82 
84 {
85  // MCParticles, MCTrackCands, MCTracks needed for this module
86  StoreArray<PXDCluster> pxdClusters;
87  pxdClusters.isRequired();
88 
89  StoreArray<SVDCluster> svdClusters;
90  svdClusters.isRequired();
91 
93  mcParticles.isRequired();
94 
96  mcTrackCands.isRequired();
97  StoreArray<genfit::TrackCand> idealMCTrackCands("idealMCTrackCands");
98  idealMCTrackCands.registerInDataStore(DataStore::c_ErrorIfAlreadyRegistered);
99 
100  StoreArray<PXDTrueHit> pxdTrueHits;
101  pxdTrueHits.isRequired();
102 
103  StoreArray<SVDTrueHit> svdTrueHits;
104  svdTrueHits.isRequired();
105 
106  //create list of histograms to be saved in the rootfile
107  m_histoList = new TList;
108 
109  //set the ROOT File
110  m_rootFilePtr = new TFile(m_rootFileName.c_str(), "RECREATE");
111 
112  //histograms to produce efficiency plots
113  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
114  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()};
115  Double_t bins_lambda[10 + 1];
116  Double_t width_lambda = TMath::Pi() / 10;
117  Double_t bins_phi[14 + 1];
118  Double_t width_phi = 2 * TMath::Pi() / 14;
119  for (int bin = 0; bin < 14 + 1; bin++)
120  bins_phi[bin] = - TMath::Pi() + bin * width_phi;
121 
122  for (int bin = 0; bin < 10 + 1; bin++) {
123  bins_lambda[bin] = - TMath::Pi() / 2 + bin * width_lambda;
124  B2DEBUG(1, bins_lambda[bin] << " " << bins_theta[bin]);
125  }
126 
127  m_h3_MCParticle = createHistogram3D("h3MCParticle", "entry per MCParticle",
128  9, bins_pt, "p_{t} (GeV/c)",
129  10, bins_lambda, "#lambda",
130  14, bins_phi, "#phi" /*, m_histoList*/);
131 
132  m_h3_idealMCTrackCand = (TH3F*)duplicateHistogram("h3idealMCTrackCand",
133  "entry per idealMCTrackCand",
134  m_h3_MCParticle /*, m_histoList*/);
135 
136  m_h3_MCTrackCand = (TH3F*)duplicateHistogram("h3MCTrackCand",
137  "entry per MCTrackCand",
138  m_h3_MCParticle /*, m_histoList*/);
139 
140 
141  m_h1_thetaMS_SVD = new TH1F("h1thetaMS_SVD", "Multiple Scattering Angle (SVD)", 500, 0, 500);
143  m_h1_thetaMS_SVD->GetXaxis()->SetTitle("#theta_{MS} (mrad)");
144 
145  m_h1_thetaMS_PXD = (TH1F*) duplicateHistogram("h1thetaMS_PXD", "Multiple Scattering Angle (PXD)", m_h1_thetaMS_SVD, m_histoList);
146 
147  m_h1_dR = new TH1F("h1dR", "dR, annulus half width", 1000, 0, 5);
148  m_histoList->Add(m_h1_dR);
149  m_h1_dR->GetXaxis()->SetTitle("dR (cm)");
150 
151  m_h1_dRoverR = new TH1F("h1dRoverR", "dR over helix radius", 1000, 0, 0.1);
153  m_h1_dRoverR->GetXaxis()->SetTitle("dR/R");
154 
155  m_h1_distOVERdR = new TH1F("h1distOVERdR", "(hit radius - helix radius)/dR", 100, -5, 5);
157  m_h1_distOVERdR->GetXaxis()->SetTitle("(hit R - helix R)/dR");
158 
159  m_h1_hitRadius_accepted = new TH1F("h1hitRadAccep", "hit radius for accepted hits", 100, 0, 500);
161  m_h1_hitRadius_accepted->GetXaxis()->SetTitle("hit radius (cm)");
162 
163  m_h1_hitRadius_rejected = new TH1F("h1hitRadRejec", "hit radius for rejected hits", 100, 0, 500);
165  m_h1_hitRadius_rejected->GetXaxis()->SetTitle("hit radius (cm)");
166 
167  m_h1_hitDistance_accepted = new TH1F("h1hitDistCentAccep", "hit distance from 0,0 for accepted hits", 100, 0, 15);
169  m_h1_hitDistance_accepted->GetXaxis()->SetTitle("hit distance (cm)");
170 
171  m_h1_hitDistance_rejected = new TH1F("h1hitDistCentRejec", "hit distance from 0,0 for rejected hits", 100, 0, 15);
173  m_h1_hitDistance_rejected->GetXaxis()->SetTitle("hit distance (cm)");
174 
175  m_h1_lapTime = new TH1F("h1LapTime", "lap time", 200, 0, 100);
177  m_h1_lapTime->GetXaxis()->SetTitle("time (ns)");
178 
179  m_h1_timeDifference = (TH1F*)duplicateHistogram("h1TimeDiff", "Hit Time Difference",
181  m_h1_diffOVERlap = new TH1F("h1HitDiffOVERlap", "Hit Time Difference over Lap Time",
182  100, 0, 1.5);
184  m_h1_diffOVERlap->GetXaxis()->SetTitle("hit difference/lap");
185 
186  m_h1_nGoodTrueHits = new TH1F("h1nTrueHitsGoods", "Number of True Hits for Accepted Tracks", 20, 0, 20);
188  m_h1_nGoodTrueHits->GetXaxis()->SetTitle("number of hits");
189 
190  m_h1_nBadTrueHits = new TH1F("h1nTrueHitsBads", "Number of True Hits for Rejected Tracks", 10, 0, 10);
192  m_h1_nBadTrueHits->GetXaxis()->SetTitle("number of hits");
193 
194  m_h1_nGood1dInfo = new TH1F("h1nGood1Dinfo", "Number of 1D Info for Accepted Tracks", 20, 0, 20);
196  m_h1_nGood1dInfo->GetXaxis()->SetTitle("number of hits");
197 
198  m_h1_nBad1dInfo = new TH1F("h1nBad1Dinfo", "Number of 1D Info for Rejected Tracks", 20, 0, 20);
200  m_h1_nBad1dInfo->GetXaxis()->SetTitle("number of hits");
201 
202  m_h1_firstRejectedHit = new TH1F("h1idealMCTCnHit", "idealMCTrackCands number of hits", 40, 0, 40);
204  m_h1_firstRejectedHit->GetXaxis()->SetTitle("# idealMCTrackCands hits");
205 
206  m_h1_firstRejectedOVERMCHit = new TH1F("h1FirstRejOVERmc", "# idealMCTrackCands hits / # MCTrackCands hits", 100, 0, 1);
208  m_h1_firstRejectedOVERMCHit->GetXaxis()->SetTitle("# idealMCTrackCands hits / # MCTrackCands hits");
209 
210  m_h1_MCTrackCandNhits = (TH1F*)duplicateHistogram("h1MCTrackCandNhits", "number of MCTrackCands hits", m_h1_firstRejectedHit,
211  m_histoList);
212 }
213 
214 
216 {
217  nWedge = 0;
218  nBarrel = 0;
219 }
220 
221 
223 {
224  B2Vector3D magField = BFieldManager::getField(0, 0, 0) / Unit::T;
225 
226  B2DEBUG(1, "+++++ 1. loop on MCTrackCands");
227 
228  StoreArray<genfit::TrackCand> idealMCTrackCands("idealMCTrackCands");
229  StoreArray<genfit::TrackCand> mcTrackCands;
230  StoreArray<PXDCluster> pxdClusters;
231  StoreArray<SVDCluster> svdClusters;
232 
234 
235  //1.a retrieve the MCTrackCands
236  BOOST_FOREACH(genfit::TrackCand & mcTrackCand, mcTrackCands) {
237 
238  int nGoodTrueHits = 0;
239  int nBadTrueHits = 0;
240  int nGood1Dinfo = 0;
241 
242  B2DEBUG(1, " a NEW MCTrackCand ");
243 
244  //1.b retrieve the MCParticle
245  RelationVector<MCParticle> MCParticles_fromMCTrackCand = DataStore::getRelationsWithObj<MCParticle>(&mcTrackCand);
246 
247  B2DEBUG(1, "~~~ " << MCParticles_fromMCTrackCand.size() << " MCParticles related to this MCTrackCand");
248  for (int mcp = 0; mcp < (int)MCParticles_fromMCTrackCand.size(); mcp++) { //should be ONE
249 
250  MCParticle mcParticle = *MCParticles_fromMCTrackCand[mcp];
251 
252  B2DEBUG(1, " a NEW charged MC Particle, " << mcParticle.getIndex() << ", " << mcParticle.getPDG());
253 
254  MCParticleInfo mcParticleInfo(mcParticle, magField);
255 
256  TVector3 decayVertex = mcParticle.getProductionVertex();
257  TVector3 mom = mcParticle.getMomentum();
258  double charge = mcParticle.getCharge();
259  double omega = mcParticleInfo.getOmega();
260  double px = mom.Px();
261  double py = mom.Py();
262  double pt = mom.Pt();
263  double x = decayVertex.X();
264  double y = decayVertex.Y();
265  double R = 1 / abs(omega); //cm
266 
267  m_h3_MCTrackCand->Fill(mcParticleInfo.getPt(), mcParticleInfo.getLambda(), mcParticleInfo.getPphi());
268 
269  double alpha = R / pt * charge; //cm/GeV
270  double Cx = x + alpha * py; //cm
271  double Cy = y - alpha * px; //cm
272 
273  TVector3 center(Cx, Cy, 0);
274 
275  //recover Clusters and loop on them
276  int Nhits = mcTrackCand.getNHits();
277  m_h1_MCTrackCandNhits->Fill(Nhits);
278 
279  int cluster = 0;
280  bool hasTrueHit = true;
281  bool isAccepted = true;
282  int firstRejectedHit = Nhits + 1;
283  double prevHitRadius = abs(1 / omega);
284 
285  double lapTime = 2 * M_PI * mcParticle.getEnergy() / 0.299792 / magField.Z();
286  double FirstHitTime = -1;
287  double HitTime = -1;
288 
289  bool isFirstSVDhit = true;
290 
291  while (cluster < Nhits && isAccepted && hasTrueHit) {
292  int detId, hitId;
293  mcTrackCand.getHit(cluster, detId, hitId);
294 
295  bool hasPXDCluster = false;
296  bool hasSVDuCluster = false;
297  bool hasSVDvCluster = false;
298 
299  double uCoor = 0;
300  double vCoor = 0;
301  VxdID sensor = 0;
302 
303  double thetaMS = 0;
304 
305  if (detId == Const::PXD && m_usePXD) {
306 
307  PXDCluster* aPXDCluster = pxdClusters[hitId];
308  RelationVector<PXDTrueHit> PXDTrueHit_fromPXDCluster = aPXDCluster->getRelationsWith<PXDTrueHit>();
309  if (PXDTrueHit_fromPXDCluster.size() == 0) {
310  B2WARNING("What's happening?!? no True Hit associated to the PXD Cluster");
311  hasTrueHit = false;
312  isAccepted = false;
313  continue;
314  }
315 
316  PXDTrueHit* aPXDTrueHit = PXDTrueHit_fromPXDCluster[0];
317  thetaMS = compute_thetaMS(mcParticleInfo, aPXDTrueHit);
318  m_h1_thetaMS_PXD->Fill(thetaMS / 2 * 1000); //PXD
319 
320  uCoor = aPXDTrueHit->getU();
321  vCoor = aPXDTrueHit->getV();
322  sensor = aPXDTrueHit->getSensorID();
323  if (cluster == 0) {
324  FirstHitTime = aPXDTrueHit->getGlobalTime();
325  HitTime = FirstHitTime;
326  } else
327  HitTime = aPXDTrueHit->getGlobalTime();
328 
329  hasPXDCluster = true;
330  } else if (detId == Const::SVD) {
331  SVDCluster* aSVDCluster = svdClusters[hitId];
332  RelationVector<SVDTrueHit> SVDTrueHit_fromSVDCluster = aSVDCluster->getRelationsWith<SVDTrueHit>();
333  if (SVDTrueHit_fromSVDCluster.size() == 0) {
334  B2WARNING("What's happening?!? no True Hit associated to the SVD Cluster");
335  hasTrueHit = false;
336  isAccepted = false;
337  continue;
338  }
339 
340  SVDTrueHit* aSVDTrueHit = SVDTrueHit_fromSVDCluster[0];
341 
342  thetaMS = compute_thetaMS(mcParticleInfo, aSVDTrueHit);
343  m_h1_thetaMS_SVD->Fill(thetaMS * 1000); //SVD
344 
345  uCoor = aSVDTrueHit->getU();
346  vCoor = aSVDTrueHit->getV();
347  sensor = aSVDTrueHit->getSensorID();
348  if (isFirstSVDhit) {
349  FirstHitTime = aSVDTrueHit->getGlobalTime();
350  HitTime = FirstHitTime;
351  isFirstSVDhit = false;
352 
353  } else
354  HitTime = aSVDTrueHit->getGlobalTime();
355  if (aSVDCluster->isUCluster())
356  hasSVDuCluster = true;
357  else
358  hasSVDvCluster = true;
359  } else {
360  cluster++;
361  continue;
362  }
363 
364  const VXD::SensorInfoBase& aSensorInfo = aGeometry.getSensorInfo(sensor);
365  bool accepted4 = true;
366  if (m_applyWedge) {
367  if (aSensorInfo.getForwardWidth() != aSensorInfo.getBackwardWidth()) {
368  nWedge++;
369  accepted4 = false;
370  } else
371  nBarrel++;
372  }
373 
374  TVector3 globalHit = aSensorInfo.pointToGlobal(TVector3(uCoor, vCoor, 0), true);
375  double hitRadius = theDistance(center, globalHit);
376 
377  bool accepted1 = true;
378  if (m_applySemiplane)
379  accepted1 = isInSemiPlane(semiPlane(decayVertex, center, globalHit), omega);
380 
381  if (accepted1) {
382  B2DEBUG(1, " semiplane: ACCEPTED");
383  } else {
384  B2DEBUG(1, " semiplane: REJECTED, next track");
385  }
386 
387  double dR = compute_dR(thetaMS, theDistance(TVector3(0, 0, 0), globalHit));
388  m_h1_dR->Fill(dR);
389  m_h1_dRoverR->Fill(dR * abs(omega));
390  m_h1_distOVERdR->Fill((hitRadius - abs(1 / omega)) / dR);
391 
392 
393  bool accepted2 = true;
394  if (m_applyAnnulus)
395  accepted2 = isInAnnulus(hitRadius, prevHitRadius, dR);
396 
397  prevHitRadius = hitRadius;
398 
399  if (accepted2) {
400  B2DEBUG(1, " annulus: ACCEPTED");
401  } else {
402  B2DEBUG(1, " annulus: REJECTED, next track");
403  }
404 
405  bool accepted3 = true;
406  if (m_applyLap)
407  accepted3 = isFirstLap(FirstHitTime, HitTime, lapTime);
408 
409  if (accepted3) {
410  B2DEBUG(1, " lapTime: ACCEPTED");
411  } else {
412  B2DEBUG(1, " lapTime: REJECTED, next track");
413  }
414 
415  if (accepted2 && accepted1 && accepted3 && accepted4) {
416  nGoodTrueHits ++;
417  m_h1_hitDistance_accepted->Fill(theDistance(TVector3(0, 0, 0), globalHit));
418  m_h1_hitRadius_accepted->Fill(hitRadius);
419  } else {
420  nBadTrueHits ++;
421  m_h1_hitDistance_rejected->Fill(theDistance(TVector3(0, 0, 0), globalHit));
422  m_h1_hitRadius_rejected->Fill(hitRadius);
423  if (m_removeBadHits)
424  firstRejectedHit = cluster;
425  isAccepted = false;
426  continue;
427  }
428 
429  if (hasPXDCluster)
430  nGood1Dinfo = +2;
431  else {
432  if (hasSVDuCluster)
433  nGood1Dinfo++;
434  if (hasSVDvCluster)
435  nGood1Dinfo++;
436  }
437  if (hasPXDCluster || hasSVDuCluster || hasSVDvCluster)
438  B2DEBUG(1, "cluster: ACCEPTED (" << nGood1Dinfo << ")");
439 
440  cluster++;
441  }//close loop on clusters
442 
443  if (nGood1Dinfo >= m_minHit) {
444  B2DEBUG(1, " idealMCTrackCand FOUND!! " << nGood1Dinfo << " 1D infos (" << nGoodTrueHits << " good true hits)");
445  m_h3_idealMCTrackCand->Fill(mcParticleInfo.getPt(), mcParticleInfo.getLambda(), mcParticleInfo.getPphi());
446  m_h1_nGoodTrueHits->Fill(nGoodTrueHits);
447  m_h1_nGood1dInfo->Fill(nGood1Dinfo);
448 
449  genfit::TrackCand* tmpTrackCand = new genfit::TrackCand(mcTrackCand);
450 
451  if ((int)firstRejectedHit <= (int)mcTrackCand.getNHits()) {
452  tmpTrackCand->reset();
453  for (int hit = 0; hit < firstRejectedHit; hit++)
454  if (mcTrackCand.getHit(hit))
455  tmpTrackCand->addHit(mcTrackCand.getHit(hit));
456  tmpTrackCand->sortHits();
457  }
458  idealMCTrackCands.appendNew(*tmpTrackCand);
459 
460  m_h1_firstRejectedHit->Fill(tmpTrackCand->getNHits());
461  m_h1_firstRejectedOVERMCHit->Fill((float)tmpTrackCand->getNHits() / mcTrackCand.getNHits());
462  } else {
463  B2DEBUG(1, " too few good hits (" << nGood1Dinfo << ") to track this one ( vs " << nGoodTrueHits << " true hits)");
464  m_h1_nBadTrueHits->Fill(nGoodTrueHits);
465  m_h1_nBad1dInfo->Fill(nGood1Dinfo);
466  }
467 
468  B2DEBUG(1, "");
469  }//close loop on MCParticles
470  }//close loop on MCTrackCands
471 }
472 
473 
475 {
476  B2INFO("** MCTrackCandClassifier parameters **");
477  B2INFO("rootfilename = " << m_rootFileName);
478  B2INFO("use PXD informations = " << m_usePXD);
479  B2INFO("--> classification criteria:");
480  if (m_applyAnnulus)
481  B2INFO(" -) |d - R| < " << m_nSigma << " dL thetaMS");
482  if (m_applySemiplane)
483  B2INFO(" -) hit in the expected semiplane");
484  if (m_applyLap)
485  B2INFO(" -) HitTime < " << m_fraction << " lap time");
486  if (m_applyWedge)
487  B2INFO(" -) hit must be in the barrel part of the VXD");
488  B2INFO("");
489 
490  double num = 0;
491  double den = 0;
492 
493  num = m_h3_idealMCTrackCand->GetEntries();
494  den = m_h3_MCTrackCand->GetEntries();
495  double efficiency = num / den ;
496  double efficiencyErr = sqrt(efficiency * (1 - efficiency)) / sqrt(den);
497 
498  B2INFO("");
499  B2INFO("~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~");
500  B2INFO("~ MCTrackCandClassifier ~ SHORT SUMMARY ~");
501  B2INFO("");
502  B2INFO(" + overall:");
503  B2INFO(" fraction of ideal MCTrackCands = (" << efficiency * 100 << " +/- " << efficiencyErr * 100 << ")% ");
504  B2INFO("");
505  B2INFO(" # idealMCTrackCand = " << num);
506  B2INFO(" # MCTrackCand = " << den);
507  B2INFO("~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~");
508  B2INFO("");
509  B2DEBUG(1, " nWedge = " << nWedge);
510  B2DEBUG(1, " nBarrel = " << nBarrel);
511 }
512 
513 
515 {
518 
519  if (m_rootFilePtr != NULL) {
520  m_rootFilePtr->cd();
521 
522  TIter nextH(m_histoList);
523  TObject* obj;
524  while ((obj = nextH()))
525  obj->Write();
526 
527  m_rootFilePtr->Close();
528  }
529 }
530 
531 
532 double MCTrackCandClassifierModule::semiPlane(TVector3 vertex, TVector3 center, TVector3 hit)
533 {
534  TVector3 err = center - vertex;
535 
536  double semiPlane = err.Y() / err.X() * hit.X() + err.Y() / err.X() * vertex.x() - vertex.Y();
537 
538  B2DEBUG(1, "");
539  B2DEBUG(1, " SEMI-PLANE defined by: y + " << err.Y() / err.X() << " x + " << err.Y() / err.X()*vertex.x() - vertex.Y() << " = 0");
540  B2DEBUG(1, " with: center(" << center.X() << "," << center.Y() << ")");
541  B2DEBUG(1, " decayV(" << vertex.X() << "," << vertex.Y() << ")");
542  B2DEBUG(1, " vector(" << err.X() << "," << err.Y() << ")");
543  B2DEBUG(1, " y SLOPE = " << semiPlane << " VS y HIT = " << hit.Y());
544  B2DEBUG(1, " HIT - SLOPE = " << - semiPlane + hit.Y());
545 
546  if (vertex.X() < center.X())
547  return hit.Y() - semiPlane;
548  else
549  return semiPlane - hit.Y();
550 }
551 
552 
553 bool MCTrackCandClassifierModule::isInSemiPlane(double semiPlane, double omega)
554 {
555  if (semiPlane * omega > 0)
556  return true;
557  else
558  return false;
559 }
560 
561 
562 double MCTrackCandClassifierModule::theDistance(TVector3 center, TVector3 hit)
563 {
564  double xSquared = TMath::Power(center.X() - hit.X(), 2);
565  double ySquared = TMath::Power(center.Y() - hit.Y(), 2);
566 
567  return TMath::Sqrt(xSquared + ySquared);
568 }
569 
570 
571 bool MCTrackCandClassifierModule::isInAnnulus(double hitDistance, double R, double dR)
572 {
573  bool accepted = false;
574 
575  B2DEBUG(1, "");
576  B2DEBUG(1, " ANNULUS defined between radii: " << R - dR << " and " << R + dR);
577  B2DEBUG(1, " hit distance = " << hitDistance);
578  B2DEBUG(1, " helix radius = " << R);
579  B2DEBUG(1, " dR = " << dR);
580 
581  if ((hitDistance > R - dR) && (hitDistance < R + dR))
582  accepted = true;
583 
584  return accepted;
585 }
586 
587 
588 bool MCTrackCandClassifierModule::isFirstLap(double FirstHitTime, double HitTime, double LapTime)
589 {
590  bool accepted = false;
591 
592  B2DEBUG(1, "");
593  B2DEBUG(1, " lapTime: " << LapTime);
594  B2DEBUG(1, " FirstHitTime = " << FirstHitTime);
595  B2DEBUG(1, " HitTime = " << HitTime);
596  B2DEBUG(1, " difference = " << HitTime - FirstHitTime);
597 
598  m_h1_lapTime->Fill(LapTime);
599  m_h1_timeDifference->Fill(HitTime - FirstHitTime);
600  m_h1_diffOVERlap->Fill((HitTime - FirstHitTime) / LapTime);
601 
602  if (HitTime - FirstHitTime < m_fraction * LapTime)
603  accepted = true;
604 
605  return accepted;
606 }
607 
608 
609 TH3F* MCTrackCandClassifierModule::createHistogram3D(const char* name, const char* title,
610  Int_t nbinsX, Double_t minX, Double_t maxX,
611  const char* titleX,
612  Int_t nbinsY, Double_t minY, Double_t maxY,
613  const char* titleY,
614  Int_t nbinsZ, Double_t minZ, Double_t maxZ,
615  const char* titleZ,
616  TList* histoList)
617 {
618  TH3F* h = new TH3F(name, title, nbinsX, minX, maxX, nbinsY, minY, maxY, nbinsZ, minZ, maxZ);
619 
620  h->GetXaxis()->SetTitle(titleX);
621  h->GetYaxis()->SetTitle(titleY);
622  h->GetZaxis()->SetTitle(titleZ);
623 
624  if (histoList)
625  histoList->Add(h);
626 
627  return h;
628 }
629 
630 
631 TH3F* MCTrackCandClassifierModule::createHistogram3D(const char* name, const char* title,
632  Int_t nbinsX, Double_t* binsX,
633  const char* titleX,
634  Int_t nbinsY, Double_t* binsY,
635  const char* titleY,
636  Int_t nbinsZ, Double_t* binsZ,
637  const char* titleZ,
638  TList* histoList)
639 {
640  TH3F* h = new TH3F(name, title, nbinsX, binsX, nbinsY, binsY, nbinsZ, binsZ);
641 
642  h->GetXaxis()->SetTitle(titleX);
643  h->GetYaxis()->SetTitle(titleY);
644  h->GetZaxis()->SetTitle(titleZ);
645 
646  if (histoList)
647  histoList->Add(h);
648 
649  return h;
650 }
651 
652 
653 TH1* MCTrackCandClassifierModule::duplicateHistogram(const char* newname, const char* newtitle,
654  TH1* h, TList* histoList)
655 {
656  TH1F* h1 = dynamic_cast<TH1F*>(h);
657  TH2F* h2 = dynamic_cast<TH2F*>(h);
658  TH3F* h3 = dynamic_cast<TH3F*>(h);
659 
660  TH1* newh = 0;
661 
662  if (h1)
663  newh = new TH1F(*h1);
664  if (h2)
665  newh = new TH2F(*h2);
666  if (h3)
667  newh = new TH3F(*h3);
668 
669  newh->SetName(newname);
670  newh->SetTitle(newtitle);
671 
672  if (histoList)
673  histoList->Add(newh);
674 
675  return newh;
676 }
677 
678 
680 {
681  //normalized to MCTrackCands
682  TH1F* h_effMCTC_pt = createHistogramsRatio("heffMCTCpt", "fraction of idealMCTrackCand VS pt", m_h3_idealMCTrackCand,
683  m_h3_MCTrackCand, true, 0);
684  histoList->Add(h_effMCTC_pt);
685 
686  TH1F* h_effMCTC_theta = createHistogramsRatio("heffMCTCtheta", "fraction of idealMCTrackCandVS #lambda", m_h3_idealMCTrackCand,
687  m_h3_MCTrackCand, true, 1);
688  histoList->Add(h_effMCTC_theta);
689 
690  TH1F* h_effMCTC_phi = createHistogramsRatio("heffMCTCphi", "fraction of idealMCTrackCand VS #phi", m_h3_idealMCTrackCand,
691  m_h3_MCTrackCand, true, 2);
692  histoList->Add(h_effMCTC_phi);
693 }
694 
695 
697 {
698  //normalized to MCTrackCands
699  TH1F* h_ineffMCTC_pt = createHistogramsRatio("hineffMCTCpt", "1 - fraction of idealMCTrackCand VS pt", m_h3_idealMCTrackCand,
700  m_h3_MCTrackCand, false, 0);
701  histoList->Add(h_ineffMCTC_pt);
702 
703  TH1F* h_ineffMCTC_theta = createHistogramsRatio("hineffMCTCtheta", "1 - fraction of idealMCTrackCandVS #lambda",
705  histoList->Add(h_ineffMCTC_theta);
706 
707  TH1F* h_ineffMCTC_phi = createHistogramsRatio("hineffMCTCphi", "1 - fraction of idealMCTrackCand VS #phi", m_h3_idealMCTrackCand,
708  m_h3_MCTrackCand, false, 2);
709  histoList->Add(h_ineffMCTC_phi);
710 }
711 
712 
713 TH1F* MCTrackCandClassifierModule::createHistogramsRatio(const char* name, const char* title,
714  TH1* hNum, TH1* hDen, bool isEffPlot,
715  int axisRef)
716 {
717  TH1F* h1den = dynamic_cast<TH1F*>(hDen);
718  TH1F* h1num = dynamic_cast<TH1F*>(hNum);
719  TH2F* h2den = dynamic_cast<TH2F*>(hDen);
720  TH2F* h2num = dynamic_cast<TH2F*>(hNum);
721  TH3F* h3den = dynamic_cast<TH3F*>(hDen);
722  TH3F* h3num = dynamic_cast<TH3F*>(hNum);
723 
724  TH1* hden = 0;
725  TH1* hnum = 0;
726 
727  if (h1den) {
728  hden = new TH1F(*h1den);
729  hnum = new TH1F(*h1num);
730  }
731  if (h2den) {
732  hden = new TH2F(*h2den);
733  hnum = new TH2F(*h2num);
734  }
735  if (h3den) {
736  hden = new TH3F(*h3den);
737  hnum = new TH3F(*h3num);
738  }
739 
740  TAxis* the_axis;
741  TAxis* the_other1;
742  TAxis* the_other2;
743 
744  if (axisRef == 0) {
745  the_axis = hden->GetXaxis();
746  the_other1 = hden->GetYaxis();
747  the_other2 = hden->GetZaxis();
748  } else if (axisRef == 1) {
749  the_axis = hden->GetYaxis();
750  the_other1 = hden->GetXaxis();
751  the_other2 = hden->GetZaxis();
752  } else if (axisRef == 2) {
753  the_axis = hden->GetZaxis();
754  the_other1 = hden->GetXaxis();
755  the_other2 = hden->GetYaxis();
756  } else
757  return NULL;
758 
759 
760  TH1F* h;
761  if (the_axis->GetXbins()->GetSize())
762  h = new TH1F(name, title, the_axis->GetNbins(), (the_axis->GetXbins())->GetArray());
763  else
764  h = new TH1F(name, title, the_axis->GetNbins(), the_axis->GetXmin(), the_axis->GetXmax());
765  h->GetXaxis()->SetTitle(the_axis->GetTitle());
766 
767  h->GetYaxis()->SetRangeUser(0.00001, 1);
768 
769  Int_t bin = 0;
770  Int_t nBins = 0;
771 
772  for (int the_bin = 1; the_bin < the_axis->GetNbins() + 1; the_bin++) {
773 
774  double num = 0;
775  double den = 0 ;
776 
777  for (int other1_bin = 1; other1_bin < the_other1->GetNbins() + 1; other1_bin++)
778  for (int other2_bin = 1; other2_bin < the_other2->GetNbins() + 1; other2_bin++) {
779 
780  if (axisRef == 0) bin = hden->GetBin(the_bin, other1_bin, other2_bin);
781  else if (axisRef == 1) bin = hden->GetBin(other1_bin, the_bin, other2_bin);
782  else if (axisRef == 2) bin = hden->GetBin(other1_bin, other2_bin, the_bin);
783 
784  if (hden->IsBinUnderflow(bin))
785  B2DEBUG(1, " bin = " << bin << "(" << the_bin << "," << other1_bin << "," << other2_bin << "), UNDERFLOW");
786  if (hden->IsBinOverflow(bin))
787  B2DEBUG(1, " bin = " << bin << "(" << the_bin << "," << other1_bin << "," << other2_bin << "), OVERFLOW");
788 
789  num += hnum->GetBinContent(bin);
790  den += hden->GetBinContent(bin);
791 
792  nBins++;
793 
794  }
795  double eff = 0;
796  double err = 0;
797 
798  if (den > 0) {
799  eff = (double)num / den;
800  err = sqrt(eff * (1 - eff)) / sqrt(den);
801  }
802 
803  if (isEffPlot) {
804  h->SetBinContent(the_bin, eff);
805  h->SetBinError(the_bin, err);
806  } else {
807  h->SetBinContent(the_bin, 1 - eff);
808  h->SetBinError(the_bin, err);
809  }
810 
811  }
812 
813  return h;
814 }
815 
816 
817 float MCTrackCandClassifierModule::compute_dR(double thetaMS, double hitDistance)
818 {
819  double dL;
820  if (hitDistance < 1.8) //L1
821  dL = 0.4;
822  else if (hitDistance < 3) //L2
823  dL = 0.8;
824  else if (hitDistance < 5) //L3
825  dL = 1.6;
826  else if (hitDistance < 9) //L4
827  dL = 4.2;
828  else if (hitDistance < 12) //L5
829  dL = 2.4;
830  else dL = 3.1;
831 
832  if ((hitDistance < 3) && (hitDistance > 1.2))
833  thetaMS = thetaMS / 2;
834 
835  double dR = m_nSigma * dL * thetaMS;
836 
837  return dR;
838 };
839 
840 
842 {
843  // double thetaMS = 0.0136 * 14 * sqrt(0.008); //SVD, PXD is half of it
844  double thetaMS = 0.0136 * 14; //SVD, PXD is half of it
845 
846  double p = mcParticleInfo.getP();
847  // double pt = mcParticleInfo.getPt();
848  double E = mcParticleInfo.getEnergy();
849 
850  double X = sqrt(pow(aTrueHit->getEntryU() - aTrueHit->getExitU(), 2) +
851  pow(aTrueHit->getEntryV() - aTrueHit->getExitV(), 2) +
852  pow(aTrueHit->getEntryW() - aTrueHit->getExitW(), 2));
853 
854  double X0 = 21.82; // g cm-2
855  double rho = 2.329; // g cm-3
856  thetaMS = thetaMS / (p * p / E) * sqrt(X / X0 * rho);
857 
858  return thetaMS;
859 };
Belle2::StoreArray::appendNew
T * appendNew()
Construct a new T object at the end of the array.
Definition: StoreArray.h:256
Belle2::RelationVector::size
size_t size() const
Get number of relations.
Definition: RelationVector.h:98
Belle2::VXD::SensorInfoBase::getForwardWidth
double getForwardWidth() const
Convinience Wrapper to return width at forward side.
Definition: SensorInfoBase.h:108
Belle2::MCTrackCandClassifierModule::m_h1_nBadTrueHits
TH1F * m_h1_nBadTrueHits
Histogram.
Definition: MCTrackCandClassifierModule.h:131
Belle2::MCParticle::getIndex
int getIndex() const
Get 1-based index of the particle in the corresponding MCParticle list.
Definition: MCParticle.h:241
Belle2::MCTrackCandClassifierModule::beginRun
void beginRun() override
Begin run.
Definition: MCTrackCandClassifierModule.cc:215
Belle2::MCParticle::getEnergy
float getEnergy() const
Return particle energy in GeV.
Definition: MCParticle.h:158
Belle2::MCTrackCandClassifierModule::isInAnnulus
bool isInAnnulus(double hitDistance, double R, double dR)
Function to check if hitDistance is within a given annulus.
Definition: MCTrackCandClassifierModule.cc:571
Belle2::MCTrackCandClassifierModule::addEfficiencyPlots
void addEfficiencyPlots(TList *graphList=NULL)
Function to create efficiency plots and add them to list.
Definition: MCTrackCandClassifierModule.cc:679
Belle2::MCParticleInfo
This struct is used by the TrackingPerformanceEvaluation Module to save information of reconstructed ...
Definition: MCParticleInfo.h:34
Belle2::MCParticle::getCharge
float getCharge() const
Return the particle charge defined in TDatabasePDG.
Definition: MCParticle.cc:35
Belle2::VxdID
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:43
Belle2::MCTrackCandClassifierModule::createHistogramsRatio
TH1F * createHistogramsRatio(const char *name, const char *title, TH1 *hNum, TH1 *hDen, bool isEffPlot, int axisRef)
Function to create a ratio histogram from two histograms.
Definition: MCTrackCandClassifierModule.cc:713
Belle2::MCTrackCandClassifierModule::m_rootFileName
std::string m_rootFileName
root file name
Definition: MCTrackCandClassifierModule.h:79
Belle2::MCTrackCandClassifierModule::m_h1_nGoodTrueHits
TH1F * m_h1_nGoodTrueHits
Histogram.
Definition: MCTrackCandClassifierModule.h:130
Belle2::MCTrackCandClassifierModule::compute_thetaMS
float compute_thetaMS(MCParticleInfo &mcParticleInfo, VXDTrueHit *aTrueHit)
Calculate thetaMS.
Definition: MCTrackCandClassifierModule.cc:841
Belle2::VXDTrueHit::getExitU
float getExitU() const
Return local u coordinate of hit at the endpoint of the track.
Definition: VXDTrueHit.h:88
Belle2::PXDTrueHit
Class PXDTrueHit - Records of tracks that either enter or leave the sensitive volume.
Definition: PXDTrueHit.h:35
Belle2::RelationsInterface::getRelationsWith
RelationVector< T > getRelationsWith(const std::string &name="", const std::string &namedRelation="") const
Get the relations between this object and another store array.
Definition: RelationsObject.h:232
Belle2::Module::setDescription
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:216
Belle2::BFieldManager::getField
static void getField(const double *pos, double *field)
return the magnetic field at a given position.
Definition: BFieldManager.h:110
Belle2::MCTrackCandClassifierModule::m_applyLap
bool m_applyLap
Whether to require that the hit belongs to the first lap in the transverse plane.
Definition: MCTrackCandClassifierModule.h:89
genfit::TrackCand
Track candidate – seed values and indices.
Definition: TrackCand.h:69
Belle2::MCTrackCandClassifierModule::m_fraction
double m_fraction
Fraction of lap.
Definition: MCTrackCandClassifierModule.h:97
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::Module::c_ParallelProcessingCertified
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:82
Belle2::VXDTrueHit::getSensorID
VxdID getSensorID() const
Return the Sensor ID.
Definition: VXDTrueHit.h:72
Belle2::VXDTrueHit::getEntryV
float getEntryV() const
Return local v coordinate of the start point of the track.
Definition: VXDTrueHit.h:84
Belle2::MCTrackCandClassifierModule::m_h1_nBad1dInfo
TH1F * m_h1_nBad1dInfo
Histogram.
Definition: MCTrackCandClassifierModule.h:133
Belle2::MCParticleInfo::getP
double getP()
Getter for magnitut of momentum.
Definition: MCParticleInfo.h:50
Belle2::SVDTrueHit
Class SVDTrueHit - Records of tracks that either enter or leave the sensitive volume.
Definition: SVDTrueHit.h:35
Belle2::VXDTrueHit
Class VXDTrueHit - Records of tracks that either enter or leave the sensitive volume.
Definition: VXDTrueHit.h:38
Belle2::VXD::SensorInfoBase
Base class to provide Sensor Information for PXD and SVD.
Definition: SensorInfoBase.h:40
Belle2::VXDTrueHit::getEntryW
float getEntryW() const
Return local w coordinate of the start point of the track.
Definition: VXDTrueHit.h:86
Belle2::MCTrackCandClassifierModule::terminate
void terminate() override
Termination action.
Definition: MCTrackCandClassifierModule.cc:514
Belle2::MCTrackCandClassifierModule::m_applySemiplane
bool m_applySemiplane
Wether to require that the hit is in the expected semiplane.
Definition: MCTrackCandClassifierModule.h:87
Belle2::MCTrackCandClassifierModule::isFirstLap
bool isFirstLap(double FirstHitTime, double HitTime, double LapTime)
Function to check if a hitTime is within a given lapTime, under consideration of m_fraction and with ...
Definition: MCTrackCandClassifierModule.cc:588
Belle2::B2Vector3::Z
DataType Z() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:434
Belle2::MCTrackCandClassifierModule::m_h1_timeDifference
TH1F * m_h1_timeDifference
Histogram.
Definition: MCTrackCandClassifierModule.h:127
Belle2::MCTrackCandClassifierModule::m_h1_nGood1dInfo
TH1F * m_h1_nGood1dInfo
Histogram.
Definition: MCTrackCandClassifierModule.h:132
Belle2::MCTrackCandClassifierModule::m_h1_hitRadius_accepted
TH1F * m_h1_hitRadius_accepted
Histogram.
Definition: MCTrackCandClassifierModule.h:117
Belle2::B2Vector3< double >
Belle2::MCTrackCandClassifierModule::m_mcTrackCandsColName
std::string m_mcTrackCandsColName
TrackCand list name.
Definition: MCTrackCandClassifierModule.h:75
Belle2::MCTrackCandClassifierModule::m_applyWedge
bool m_applyWedge
Whether to require that the hit belong to the barrel part of the SVD.
Definition: MCTrackCandClassifierModule.h:91
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::MCTrackCandClassifierModule::m_rootFilePtr
TFile * m_rootFilePtr
Pointer to root file used for storing histograms.
Definition: MCTrackCandClassifierModule.h:105
Belle2::MCTrackCandClassifierModule::m_h3_MCTrackCand
TH3F * m_h3_MCTrackCand
Histogram.
Definition: MCTrackCandClassifierModule.h:110
Belle2::Module::setPropertyFlags
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:210
Belle2::Unit::T
static const double T
[tesla]
Definition: Unit.h:130
Belle2::MCParticleInfo::getPt
double getPt()
Getter for transverse momentum.
Definition: MCParticleInfo.h:48
Belle2::MCTrackCandClassifierModule::m_h1_firstRejectedOVERMCHit
TH1F * m_h1_firstRejectedOVERMCHit
Histogram.
Definition: MCTrackCandClassifierModule.h:124
Belle2::MCTrackCandClassifierModule::m_h1_diffOVERlap
TH1F * m_h1_diffOVERlap
Histogram.
Definition: MCTrackCandClassifierModule.h:128
Belle2::MCTrackCandClassifierModule::m_h1_hitDistance_accepted
TH1F * m_h1_hitDistance_accepted
Histogram.
Definition: MCTrackCandClassifierModule.h:119
Belle2::RelationVector
Class for type safe access to objects that are referred to in relations.
Definition: DataStore.h:38
Belle2::MCTrackCandClassifierModule::m_h1_dRoverR
TH1F * m_h1_dRoverR
Histogram.
Definition: MCTrackCandClassifierModule.h:115
Belle2::VXD::GeoCache::getInstance
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:215
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::MCTrackCandClassifierModule::initialize
void initialize() override
Initializes the Module.
Definition: MCTrackCandClassifierModule.cc:83
Belle2::MCTrackCandClassifierModule::m_h1_MCTrackCandNhits
TH1F * m_h1_MCTrackCandNhits
Histogram.
Definition: MCTrackCandClassifierModule.h:122
Belle2::VXDTrueHit::getV
float getV() const
Return local v coordinate of hit.
Definition: VXDTrueHit.h:78
Belle2::MCParticleInfo::getPphi
double getPphi()
Getter for phi of momentum vector.
Definition: MCParticleInfo.h:64
Belle2::MCTrackCandClassifierModule::m_nSigma
int m_nSigma
nSigma dR
Definition: MCTrackCandClassifierModule.h:81
Belle2::MCParticleInfo::getOmega
double getOmega()
Getter for Omega.
Definition: MCParticleInfo.cc:69
Belle2::VXDTrueHit::getExitV
float getExitV() const
Return local v coordinate of hit at the endpoint of the track.
Definition: VXDTrueHit.h:90
Belle2::MCParticle::getProductionVertex
TVector3 getProductionVertex() const
Return production vertex position.
Definition: MCParticle.h:200
Belle2::MCParticle::getPDG
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:123
Belle2::MCTrackCandClassifierModule::compute_dR
float compute_dR(double thetaMS, double omega)
Calculate dR.
Definition: MCTrackCandClassifierModule.cc:817
Belle2::MCTrackCandClassifierModule::isInSemiPlane
bool isInSemiPlane(double semiPlane, double omega)
Function to check if a omega value is in a given semiPlane.
Definition: MCTrackCandClassifierModule.cc:553
Belle2::MCTrackCandClassifierModule::semiPlane
double semiPlane(TVector3 vertex, TVector3 center, TVector3 hit)
Function to get semiplane.
Definition: MCTrackCandClassifierModule.cc:532
Belle2::MCTrackCandClassifierModule::m_removeBadHits
bool m_removeBadHits
Whether to remove the clusters that do not satisfy the criteria from the idealMCTrackCands.
Definition: MCTrackCandClassifierModule.h:93
genfit::TrackCand::sortHits
void sortHits()
Sort the hits that were already added to the trackCand using the sorting parameters.
Definition: TrackCand.cc:222
Belle2::MCTrackCandClassifierModule::m_h1_hitDistance_rejected
TH1F * m_h1_hitDistance_rejected
Histogram.
Definition: MCTrackCandClassifierModule.h:120
Belle2::PXDCluster
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
Definition: PXDCluster.h:41
Belle2::MCTrackCandClassifierModule::m_h1_distOVERdR
TH1F * m_h1_distOVERdR
Histogram.
Definition: MCTrackCandClassifierModule.h:116
Belle2::MCTrackCandClassifierModule::m_usePXD
bool m_usePXD
Whether to use PXD.
Definition: MCTrackCandClassifierModule.h:83
Belle2::VXD::SensorInfoBase::pointToGlobal
TVector3 pointToGlobal(const TVector3 &local, bool reco=false) const
Convert a point from local to global coordinates.
Definition: SensorInfoBase.h:373
Belle2::SVDCluster::isUCluster
bool isUCluster() const
Get the direction of strips.
Definition: SVDCluster.h:118
Belle2::VXDTrueHit::getEntryU
float getEntryU() const
Return local u coordinate of hit when entering silicon.
Definition: VXDTrueHit.h:82
Belle2::MCTrackCandClassifierModule::m_h1_dR
TH1F * m_h1_dR
Histogram.
Definition: MCTrackCandClassifierModule.h:114
Belle2::MCTrackCandClassifierModule::m_mcParticlesName
std::string m_mcParticlesName
MCParticle list name.
Definition: MCTrackCandClassifierModule.h:77
Belle2::MCTrackCandClassifierModule::m_h1_hitRadius_rejected
TH1F * m_h1_hitRadius_rejected
Histogram.
Definition: MCTrackCandClassifierModule.h:118
Belle2::MCTrackCandClassifierModule::m_applyAnnulus
bool m_applyAnnulus
Whether to require that the hit is in the expected annulus.
Definition: MCTrackCandClassifierModule.h:85
Belle2::MCTrackCandClassifierModule::m_minHit
int m_minHit
Minimum number of 1D Clusters to classify the MCTrackCand as ideal.
Definition: MCTrackCandClassifierModule.h:95
Belle2::MCTrackCandClassifierModule::theDistance
double theDistance(TVector3 center, TVector3 hit)
Get distance between two points.
Definition: MCTrackCandClassifierModule.cc:562
Belle2::DataStore::c_ErrorIfAlreadyRegistered
@ c_ErrorIfAlreadyRegistered
If the object/array was already registered, produce an error (aborting initialisation).
Definition: DataStore.h:74
Belle2::SVDCluster
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
Definition: SVDCluster.h:38
Belle2::MCParticleInfo::getEnergy
double getEnergy()
Getter for energy.
Definition: MCParticleInfo.h:52
Belle2::Module::addParam
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:562
Belle2::VXD::GeoCache
Class to faciliate easy access to sensor information of the VXD like coordinate transformations or pi...
Definition: GeoCache.h:41
Belle2::MCTrackCandClassifierModule::m_h1_lapTime
TH1F * m_h1_lapTime
Histogram.
Definition: MCTrackCandClassifierModule.h:126
Belle2::MCTrackCandClassifierModule::duplicateHistogram
TH1 * duplicateHistogram(const char *newname, const char *newtitle, TH1 *h, TList *histoList=NULL)
Function to clone a histogram.
Definition: MCTrackCandClassifierModule.cc:653
Belle2::MCTrackCandClassifierModule::m_h3_idealMCTrackCand
TH3F * m_h3_idealMCTrackCand
Histogram.
Definition: MCTrackCandClassifierModule.h:109
Belle2::MCTrackCandClassifierModule::createHistogram3D
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=NULL)
Create a 3D ROOT Histogram.
Definition: MCTrackCandClassifierModule.cc:609
Belle2::MCTrackCandClassifierModule::m_h1_thetaMS_SVD
TH1F * m_h1_thetaMS_SVD
Histogram.
Definition: MCTrackCandClassifierModule.h:112
Belle2::MCParticle::getMomentum
TVector3 getMomentum() const
Return momentum.
Definition: MCParticle.h:209
Belle2::VXD::GeoCache::getSensorInfo
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:68
Belle2::MCParticle
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:43
Belle2::VXDTrueHit::getU
float getU() const
Return local u coordinate of hit.
Definition: VXDTrueHit.h:76
Belle2::StoreArray< PXDCluster >
Belle2::MCTrackCandClassifierModule::m_h1_thetaMS_PXD
TH1F * m_h1_thetaMS_PXD
Histogram.
Definition: MCTrackCandClassifierModule.h:113
Belle2::VXD::SensorInfoBase::getBackwardWidth
double getBackwardWidth() const
Convinience Wrapper to return width at backward side.
Definition: SensorInfoBase.h:100
Belle2::VXDTrueHit::getExitW
float getExitW() const
Return local w coordinate of hit at the endpoint of the track.
Definition: VXDTrueHit.h:92
genfit::TrackCand::reset
void reset()
Delete and clear the TrackCandHits.
Definition: TrackCand.cc:172
Belle2::MCParticleInfo::getLambda
double getLambda()
Getter for Lambda.
Definition: MCParticleInfo.cc:90
Belle2::VXDTrueHit::getGlobalTime
float getGlobalTime() const
Return the time when the track reached its midpoint.
Definition: VXDTrueHit.h:96
Belle2::MCTrackCandClassifierModule::nBarrel
int nBarrel
Counter for hits on barrel sensors.
Definition: MCTrackCandClassifierModule.h:102
Belle2::MCTrackCandClassifierModule::addInefficiencyPlots
void addInefficiencyPlots(TList *graphList=NULL)
Function to create inefficiency plots and add them to list.
Definition: MCTrackCandClassifierModule.cc:696
Belle2::MCTrackCandClassifierModule::nWedge
int nWedge
Counter for hits on wedged sensors.
Definition: MCTrackCandClassifierModule.h:100
Belle2::MCTrackCandClassifierModule::m_histoList
TList * m_histoList
List of histograms.
Definition: MCTrackCandClassifierModule.h:107
Belle2::MCTrackCandClassifierModule::endRun
void endRun() override
End run.
Definition: MCTrackCandClassifierModule.cc:474
Belle2::MCTrackCandClassifierModule::event
void event() override
Event function.
Definition: MCTrackCandClassifierModule.cc:222
Belle2::MCTrackCandClassifierModule::m_h3_MCParticle
TH3F * m_h3_MCParticle
Histogram.
Definition: MCTrackCandClassifierModule.h:108
Belle2::MCTrackCandClassifierModule::m_h1_firstRejectedHit
TH1F * m_h1_firstRejectedHit
Histogram.
Definition: MCTrackCandClassifierModule.h:123