Belle II Software  release-05-02-19
KLMClusterEfficiencyModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2016 Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Kirill Chilikin *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 /* Own header. */
12 #include <klm/modules/KLMClusterEfficiency/KLMClusterEfficiencyModule.h>
13 
14 /* KLM headers. */
15 #include <klm/dataobjects/bklm/BKLMHit2d.h>
16 #include <klm/dataobjects/eklm/EKLMHit2d.h>
17 
18 /* ROOT headers. */
19 #include <TCanvas.h>
20 #include <TH2F.h>
21 #include <TMarker.h>
22 #include <TStyle.h>
23 
24 using namespace Belle2;
25 
26 REG_MODULE(KLMClusterEfficiency)
27 
29  m_OutputFile(nullptr), m_ClusterTree(nullptr), m_ReconstructionTree(nullptr),
30  m_DecayVertexX(0), m_DecayVertexY(0), m_DecayVertexZ(0),
31  m_MaxDecayVertexHitAngle(0), m_ClusterX(0), m_ClusterY(0), m_ClusterZ(0),
32  m_MaxClusterHitAngle(0), m_KL0Clusters(0), m_PartlyKL0Clusters(0),
33  m_OtherClusters(0), m_NonreconstructedKL0(0),
34  m_ReconstructedKL01Cluster{0, 0, 0}, m_ExactlyReconstructedKL0(0),
35  m_ReconstructedKL02Clusters{0, 0, 0, 0, 0, 0},
36  m_ReconstructedKL03Clusters(0)
37 {
38  setDescription("Module for KLM cluster reconstruction efficiency studies.");
39  addParam("SaveClusterData", m_SaveClusterData,
40  "Whether to save cluster data or not.", false);
41  addParam("SaveReconstructionData", m_SaveReconstructionData,
42  "Whether to save reconstruction data or not.", true);
43  addParam("OutputFile", m_OutputFileName, "Output file.",
44  std::string("KLMClusterEfficiency.root"));
45  addParam("EventsClusterHistograms", m_EventsClusterHistograms,
46  "Draw cluster histograms for this number of events.", 0);
47 }
48 
50 {
51 }
52 
54 {
56  m_OutputFile = new TFile(m_OutputFileName.c_str(), "recreate");
57  if (m_SaveClusterData) {
58  m_ClusterTree = new TTree("klm_cluster", "");
59  m_ClusterTree->Branch("DecayVertexX", &m_DecayVertexX, "DecayVertexX/F");
60  m_ClusterTree->Branch("DecayVertexY", &m_DecayVertexY, "DecayVertexY/F");
61  m_ClusterTree->Branch("DecayVertexZ", &m_DecayVertexZ, "DecayVertexZ/F");
62  m_ClusterTree->Branch("MaxDecayVertexHitAngle", &m_MaxDecayVertexHitAngle,
63  "MaxDecayVertexHitAngle/F");
64  m_ClusterTree->Branch("ClusterX", &m_ClusterX, "ClusterX/F");
65  m_ClusterTree->Branch("ClusterY", &m_ClusterY, "ClusterY/F");
66  m_ClusterTree->Branch("ClusterZ", &m_ClusterZ, "ClusterZ/F");
67  m_ClusterTree->Branch("MaxClusterHitAngle", &m_MaxClusterHitAngle,
68  "MaxClusterHitAngle/F");
69  }
71  m_ReconstructionTree = new TTree("klm_reconstruction", "");
72  m_ReconstructionTree->Branch("KL0Clusters", &m_KL0Clusters,
73  "KL0Clusters/I");
74  m_ReconstructionTree->Branch("PartlyKL0Clusters", &m_PartlyKL0Clusters,
75  "PartlyKL0Clusters/I");
76  m_ReconstructionTree->Branch("OtherClusters", &m_OtherClusters,
77  "OtherClusters/I");
78  m_ReconstructionTree->Branch("NonreconstructedKL0", &m_NonreconstructedKL0,
79  "NonreconstructedKL0/I");
80  m_ReconstructionTree->Branch("ReconstructedKL01Cluster",
82  "ReconstructedKL01Cluster[3]/I");
83  m_ReconstructionTree->Branch("ReconstructedKL02Clusters",
85  "ReconstructedKL02Clusters[6]/I");
86  m_ReconstructionTree->Branch("ReconstructedKL03Clusters",
88  "ReconstructedKL03Clusters/I");
89  }
90  gStyle->SetOptStat(0);
91 }
92 
94 {
95 }
96 
98 {
99  static int nevent = 0;
100  /* cppcheck-suppress variableScope */
101  char str[128];
102  int i1, i2, i3, n1, n2, n3;
103  int bs1, bs2, es1, es2;
104  TVector3 decayVertex, clusterPosition, hitPosition;
105  float angle;
106  /* cppcheck-suppress variableScope */
107  bool haveKL0;
108  n1 = m_KLMClusters.getEntries();
109  for (i1 = 0; i1 < n1; i1++) {
110  haveKL0 = false;
111  RelationVector<MCParticle> clusterMCParticles =
112  m_KLMClusters[i1]->getRelationsTo<MCParticle>();
113  n2 = clusterMCParticles.size();
114  for (i2 = 0; i2 < n2; i2++) {
115  if (clusterMCParticles[i2]->getPDG() == 130)
116  haveKL0 = true;
117  }
118  if (haveKL0) {
119  if (n2 == 1)
120  m_KL0Clusters++;
121  else
123  } else
124  m_OtherClusters++;
125  }
126  if (nevent < m_EventsClusterHistograms) {
127  static TH2F* hzx = new TH2F("hzx", "", 100, -300, 420, 100, -320, 320);
128  static TH2F* hzy = new TH2F("hzy", "", 100, -300, 420, 100, -320, 320);
129  static TH2F* hxy = new TH2F("hxy", "", 100, -320, 320, 100, -320, 320);
130  static TCanvas* c1 = new TCanvas();
131  static TMarker* clusterMarker = new TMarker(0, 0, 20);
132  static TMarker* hitMarker = new TMarker(0, 0, 21);
133  if (nevent == 0) {
134  hzx->GetXaxis()->SetTitle("z, cm");
135  hzx->GetYaxis()->SetTitle("x, cm");
136  hzy->GetXaxis()->SetTitle("z, cm");
137  hzy->GetYaxis()->SetTitle("y, cm");
138  hxy->GetXaxis()->SetTitle("x, cm");
139  hxy->GetYaxis()->SetTitle("y, cm");
140  }
141  hzx->Draw();
142  for (i1 = 0; i1 < n1; i1++) {
143  clusterMarker->SetMarkerColor(i1 + 1);
144  hitMarker->SetMarkerColor(i1 + 1);
145  clusterPosition = m_KLMClusters[i1]->getClusterPosition();
146  clusterMarker->DrawMarker(clusterPosition.Z(), clusterPosition.X());
147  RelationVector<BKLMHit2d> bklmHit2ds =
148  m_KLMClusters[i1]->getRelationsTo<BKLMHit2d>();
149  n2 = bklmHit2ds.size();
150  for (i2 = 0; i2 < n2; i2++) {
151  hitPosition = bklmHit2ds[i2]->getGlobalPosition();
152  hitMarker->DrawMarker(hitPosition.Z(), hitPosition.X());
153  }
154  RelationVector<EKLMHit2d> eklmHit2ds =
155  m_KLMClusters[i1]->getRelationsTo<EKLMHit2d>();
156  n2 = eklmHit2ds.size();
157  for (i2 = 0; i2 < n2; i2++) {
158  hitPosition = eklmHit2ds[i2]->getPosition();
159  hitMarker->DrawMarker(hitPosition.Z(), hitPosition.X());
160  }
161  }
162  snprintf(str, 128, "clusters%dzx.eps", nevent);
163  c1->Print(str);
164  hzy->Draw();
165  for (i1 = 0; i1 < n1; i1++) {
166  clusterMarker->SetMarkerColor(i1 + 1);
167  hitMarker->SetMarkerColor(i1 + 1);
168  clusterPosition = m_KLMClusters[i1]->getClusterPosition();
169  clusterMarker->DrawMarker(clusterPosition.Z(), clusterPosition.Y());
170  RelationVector<BKLMHit2d> bklmHit2ds =
171  m_KLMClusters[i1]->getRelationsTo<BKLMHit2d>();
172  n2 = bklmHit2ds.size();
173  for (i2 = 0; i2 < n2; i2++) {
174  hitPosition = bklmHit2ds[i2]->getGlobalPosition();
175  hitMarker->DrawMarker(hitPosition.Z(), hitPosition.Y());
176  }
177  RelationVector<EKLMHit2d> eklmHit2ds =
178  m_KLMClusters[i1]->getRelationsTo<EKLMHit2d>();
179  n2 = eklmHit2ds.size();
180  for (i2 = 0; i2 < n2; i2++) {
181  hitPosition = eklmHit2ds[i2]->getPosition();
182  hitMarker->DrawMarker(hitPosition.Z(), hitPosition.Y());
183  }
184  }
185  snprintf(str, 128, "clusters%dzy.eps", nevent);
186  c1->Print(str);
187  hxy->Draw();
188  for (i1 = 0; i1 < n1; i1++) {
189  clusterMarker->SetMarkerColor(i1 + 1);
190  hitMarker->SetMarkerColor(i1 + 1);
191  clusterPosition = m_KLMClusters[i1]->getClusterPosition();
192  clusterMarker->DrawMarker(clusterPosition.X(), clusterPosition.Y());
193  RelationVector<BKLMHit2d> bklmHit2ds =
194  m_KLMClusters[i1]->getRelationsTo<BKLMHit2d>();
195  n2 = bklmHit2ds.size();
196  for (i2 = 0; i2 < n2; i2++) {
197  hitPosition = bklmHit2ds[i2]->getGlobalPosition();
198  hitMarker->DrawMarker(hitPosition.X(), hitPosition.Y());
199  }
200  RelationVector<EKLMHit2d> eklmHit2ds =
201  m_KLMClusters[i1]->getRelationsTo<EKLMHit2d>();
202  n2 = eklmHit2ds.size();
203  for (i2 = 0; i2 < n2; i2++) {
204  hitPosition = eklmHit2ds[i2]->getPosition();
205  hitMarker->DrawMarker(hitPosition.X(), hitPosition.Y());
206  }
207  }
208  snprintf(str, 128, "clusters%dxy.eps", nevent);
209  c1->Print(str);
210  nevent++;
211  }
212  n1 = m_MCParticles.getEntries();
213  for (i1 = 0; i1 < n1; i1++) {
214  if (m_MCParticles[i1]->getPDG() != 130)
215  continue;
216  decayVertex = m_MCParticles[i1]->getDecayVertex();
217  m_DecayVertexX = decayVertex.X();
218  m_DecayVertexY = decayVertex.Y();
219  m_DecayVertexZ = decayVertex.Z();
221  RelationVector<BKLMHit2d> mcBKLMHit2ds =
222  m_MCParticles[i1]->getRelationsFrom<BKLMHit2d>();
223  n2 = mcBKLMHit2ds.size();
224  for (i2 = 0; i2 < n2; i2++) {
225  hitPosition = mcBKLMHit2ds[i2]->getGlobalPosition();
226  angle = decayVertex.Angle(hitPosition);
227  if (angle > m_MaxDecayVertexHitAngle)
228  m_MaxDecayVertexHitAngle = angle;
229  }
230  RelationVector<EKLMHit2d> mcEKLMHit2ds =
231  m_MCParticles[i1]->getRelationsFrom<EKLMHit2d>();
232  n2 = mcEKLMHit2ds.size();
233  for (i2 = 0; i2 < n2; i2++) {
234  hitPosition = mcEKLMHit2ds[i2]->getPosition();
235  angle = decayVertex.Angle(hitPosition);
236  if (angle > m_MaxDecayVertexHitAngle)
237  m_MaxDecayVertexHitAngle = angle;
238  }
239  RelationVector<KLMCluster> kl0Clusters =
240  m_MCParticles[i1]->getRelationsFrom<KLMCluster>();
241  n2 = kl0Clusters.size();
242  if (n2 == 0)
244  else if (n2 == 1) {
245  RelationVector<BKLMHit2d> bklmHit2ds =
246  kl0Clusters[0]->getRelationsTo<BKLMHit2d>();
247  RelationVector<EKLMHit2d> eklmHit2ds =
248  kl0Clusters[0]->getRelationsTo<EKLMHit2d>();
249  bs1 = bklmHit2ds.size();
250  es1 = eklmHit2ds.size();
251  if (bs1 > 0) {
252  if (es1 > 0)
254  else
256  } else if (es1 > 0)
258  RelationVector<MCParticle> m_MCParticles2 =
259  kl0Clusters[0]->getRelationsTo<MCParticle>();
260  if (m_MCParticles2.size() == 1) {
261  if (m_MCParticles2.weight(0) == 1)
263  }
264  } else if (n2 == 2) {
265  RelationVector<BKLMHit2d> bklmHit2ds1 =
266  kl0Clusters[0]->getRelationsTo<BKLMHit2d>();
267  RelationVector<BKLMHit2d> bklmHit2ds2 =
268  kl0Clusters[1]->getRelationsTo<BKLMHit2d>();
269  RelationVector<EKLMHit2d> eklmHit2ds1 =
270  kl0Clusters[0]->getRelationsTo<EKLMHit2d>();
271  RelationVector<EKLMHit2d> eklmHit2ds2 =
272  kl0Clusters[1]->getRelationsTo<EKLMHit2d>();
273  bs1 = bklmHit2ds1.size();
274  bs2 = bklmHit2ds2.size();
275  es1 = eklmHit2ds1.size();
276  es2 = eklmHit2ds2.size();
277  if (bs1 > 0 && bs2 > 0) {
278  if (es1 > 0 && es2 > 0) {
280  } else if ((es1 > 0 && es2 == 0) || (es1 == 0 && es2 > 0)) {
282  } else if (es1 == 0 && es2 == 0) {
284  }
285  } else if (bs1 > 0 && bs2 == 0) {
286  if (es1 > 0 && es2 > 0) {
288  } else if (es1 == 0 && es2 > 0) {
290  }
291  } else if (bs1 == 0 && bs2 > 0) {
292  if (es1 > 0 && es2 > 0) {
294  } else if (es1 > 0 && es2 == 0) {
296  }
297  } else if (bs1 == 0 && bs2 == 0) {
298  if (es1 > 0 && es2 > 0) {
300  }
301  }
302  } else if (n2 >= 3)
304  if (m_SaveClusterData) {
305  for (i2 = 0; i2 < n2; i2++) {
306  clusterPosition = kl0Clusters[i2]->getClusterPosition();
307  m_ClusterX = clusterPosition.X();
308  m_ClusterY = clusterPosition.Y();
309  m_ClusterZ = clusterPosition.Z();
311  RelationVector<BKLMHit2d> bklmHit2ds =
312  kl0Clusters[i2]->getRelationsTo<BKLMHit2d>();
313  n3 = bklmHit2ds.size();
314  for (i3 = 0; i3 < n3; i3++) {
315  hitPosition = bklmHit2ds[i3]->getGlobalPosition();
316  angle = clusterPosition.Angle(hitPosition);
317  if (angle > m_MaxClusterHitAngle)
318  m_MaxClusterHitAngle = angle;
319  }
320  RelationVector<EKLMHit2d> eklmHit2ds =
321  kl0Clusters[i2]->getRelationsTo<EKLMHit2d>();
322  n3 = eklmHit2ds.size();
323  for (i3 = 0; i3 < n3; i3++) {
324  hitPosition = eklmHit2ds[i3]->getPosition();
325  angle = clusterPosition.Angle(hitPosition);
326  if (angle > m_MaxClusterHitAngle)
327  m_MaxClusterHitAngle = angle;
328  }
329  m_ClusterTree->Fill();
330  }
331  }
332  }
333 }
334 
336 {
337 }
338 
340 {
341  int i, rec1Cluster, rec2Clusters;
342  rec1Cluster = 0;
343  rec2Clusters = 0;
344  for (i = 0; i < 3; i++)
345  rec1Cluster = rec1Cluster + m_ReconstructedKL01Cluster[i];
346  for (i = 0; i < 6; i++)
347  rec2Clusters = rec2Clusters + m_ReconstructedKL02Clusters[i];
348  /* Always printed once, not necessary to use LogVar. */
349  B2INFO("Total number of KLM clusters: " << m_KL0Clusters +
351  B2INFO("K_L0 clusters: " << m_KL0Clusters);
352  B2INFO("(K_L0+other) clusters: " << m_PartlyKL0Clusters);
353  B2INFO("Other clusters: " << m_OtherClusters);
354  B2INFO("Total number of generated K_L0: " << m_NonreconstructedKL0 +
355  rec1Cluster + rec2Clusters + m_ReconstructedKL03Clusters);
356  B2INFO("Nonreconstructed K_L0: " << m_NonreconstructedKL0);
357  B2INFO("K_L0 reconstructed as 1 cluster (total): " << rec1Cluster);
358  B2INFO("K_L0 reconstructed as 1 cluster (BKLM): " <<
360  B2INFO("K_L0 reconstructed as 1 cluster (BKLM/EKLM): " <<
362  B2INFO("K_L0 reconstructed as 1 cluster (EKLM): " <<
364  B2INFO("K_L0 reconstructed as 1 cluster (exact reconstruction): " <<
366  B2INFO("K_L0 reconstructed as 2 clusters (total): " << rec2Clusters);
367  B2INFO("K_L0 reconstructed as 2 clusters (2 * BKLM): " <<
369  B2INFO("K_L0 reconstructed as 2 clusters (BKLM + BKLM/EKLM): " <<
371  B2INFO("K_L0 reconstructed as 2 clusters (BKLM + EKLM): " <<
373  B2INFO("K_L0 reconstructed as 2 clusters (2 * BKLM/EKLM): " <<
375  B2INFO("K_L0 reconstructed as 2 clusters (BKLM/EKLM + EKLM): " <<
377  B2INFO("K_L0 reconstructed as 2 clusters (2 * EKLM): " <<
379  B2INFO("K_L0 reconstructed as 3 or more clusters: " <<
382  m_ReconstructionTree->Fill();
383  m_OutputFile->cd();
384  m_ReconstructionTree->Write();
385  delete m_ReconstructionTree;
386  }
387  if (m_SaveClusterData) {
388  m_OutputFile->cd();
389  m_ClusterTree->Write();
390  delete m_ClusterTree;
391  }
393  delete m_OutputFile;
394 }
395 
Belle2::RelationVector::size
size_t size() const
Get number of relations.
Definition: RelationVector.h:98
Belle2::KLMClusterEfficiencyModule::m_ReconstructedKL03Clusters
int m_ReconstructedKL03Clusters
Number of K_L0 reconstructed as >= 2 clusters.
Definition: KLMClusterEfficiencyModule.h:158
Belle2::KLMClusterEfficiencyModule::~KLMClusterEfficiencyModule
~KLMClusterEfficiencyModule()
Destructor.
Definition: KLMClusterEfficiencyModule.cc:49
Belle2::KLMClusterEfficiencyModule::m_ReconstructionTree
TTree * m_ReconstructionTree
Reconstruction tree.
Definition: KLMClusterEfficiencyModule.h:100
Belle2::KLMClusterEfficiencyModule::beginRun
void beginRun() override
Called when entering a new run.
Definition: KLMClusterEfficiencyModule.cc:93
Belle2::KLMClusterEfficiencyModule::m_OutputFileName
std::string m_OutputFileName
Output file name.
Definition: KLMClusterEfficiencyModule.h:88
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::KLMClusterEfficiencyModule::m_DecayVertexX
float m_DecayVertexX
MCParticle decay vertex X coordinate.
Definition: KLMClusterEfficiencyModule.h:103
Belle2::KLMClusterEfficiencyModule::m_MaxClusterHitAngle
float m_MaxClusterHitAngle
Maximal angle between KLM cluster and its hits.
Definition: KLMClusterEfficiencyModule.h:124
Belle2::KLMClusterEfficiencyModule::m_ExactlyReconstructedKL0
int m_ExactlyReconstructedKL0
Number of K_L0 reconstructed as 1 cluster, and this cluster should be related only to 1 (K_L0 ) MCPar...
Definition: KLMClusterEfficiencyModule.h:148
Belle2::KLMClusterEfficiencyModule::m_OtherClusters
int m_OtherClusters
Number of clusters from other particles.
Definition: KLMClusterEfficiencyModule.h:133
Belle2::KLMClusterEfficiencyModule::m_SaveReconstructionData
bool m_SaveReconstructionData
Whether to save reconstruction data or not.
Definition: KLMClusterEfficiencyModule.h:85
Belle2::KLMClusterEfficiencyModule::m_KL0Clusters
int m_KL0Clusters
Number of clusters from a K_L0.
Definition: KLMClusterEfficiencyModule.h:127
Belle2::KLMClusterEfficiencyModule::m_ClusterZ
float m_ClusterZ
Cluster Z coordinate.
Definition: KLMClusterEfficiencyModule.h:121
Belle2::KLMClusterEfficiencyModule::m_KLMClusters
StoreArray< KLMCluster > m_KLMClusters
KLM clusters.
Definition: KLMClusterEfficiencyModule.h:161
Belle2::KLMClusterEfficiencyModule::m_NonreconstructedKL0
int m_NonreconstructedKL0
Number of nonreconstructed K_L0.
Definition: KLMClusterEfficiencyModule.h:136
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::KLMClusterEfficiencyModule::m_DecayVertexY
float m_DecayVertexY
MCParticle decay vertex Y coordinate.
Definition: KLMClusterEfficiencyModule.h:106
Belle2::RelationVector
Class for type safe access to objects that are referred to in relations.
Definition: DataStore.h:38
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::KLMClusterEfficiencyModule::m_SaveClusterData
bool m_SaveClusterData
Whether to save cluster data or not.
Definition: KLMClusterEfficiencyModule.h:82
Belle2::KLMClusterEfficiencyModule::m_PartlyKL0Clusters
int m_PartlyKL0Clusters
Number of clusters from a K_L0 + other particles.
Definition: KLMClusterEfficiencyModule.h:130
Belle2::KLMClusterEfficiencyModule::m_MaxDecayVertexHitAngle
float m_MaxDecayVertexHitAngle
Maximal angle between MCParticle decay vertex and its hits.
Definition: KLMClusterEfficiencyModule.h:112
Belle2::KLMClusterEfficiencyModule::m_EventsClusterHistograms
int m_EventsClusterHistograms
Draw cluster histograms for this number of events.
Definition: KLMClusterEfficiencyModule.h:91
Belle2::KLMClusterEfficiencyModule::m_OutputFile
TFile * m_OutputFile
Output file.
Definition: KLMClusterEfficiencyModule.h:94
Belle2::KLMCluster
KLM cluster data.
Definition: KLMCluster.h:38
Belle2::KLMClusterEfficiencyModule::initialize
void initialize() override
Initializer.
Definition: KLMClusterEfficiencyModule.cc:53
Belle2::RelationVector::weight
float weight(int index) const
Get weight with index.
Definition: RelationVector.h:120
Belle2::KLMClusterEfficiencyModule::m_ClusterX
float m_ClusterX
Cluster X coordinate.
Definition: KLMClusterEfficiencyModule.h:115
Belle2::KLMClusterEfficiencyModule::m_ClusterTree
TTree * m_ClusterTree
Cluster data tree.
Definition: KLMClusterEfficiencyModule.h:97
Belle2::KLMClusterEfficiencyModule::m_ClusterY
float m_ClusterY
Cluster Y coordinate.
Definition: KLMClusterEfficiencyModule.h:118
Belle2::MCParticle
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:43
Belle2::KLMClusterEfficiencyModule::m_ReconstructedKL02Clusters
int m_ReconstructedKL02Clusters[6]
Number of K_L0 reconstructed as 2 clusters in (2 * BKLM, BKLM + BKLM/EKLM, BKLM + EKLM,...
Definition: KLMClusterEfficiencyModule.h:155
Belle2::KLMClusterEfficiencyModule::m_ReconstructedKL01Cluster
int m_ReconstructedKL01Cluster[3]
Number of K_L0 reconstructed as 1 cluster in (BKLM, BKLM/EKLM, EKLM).
Definition: KLMClusterEfficiencyModule.h:142
Belle2::KLMClusterEfficiencyModule
Module for KLM cluster reconstruction efficiency studies.
Definition: KLMClusterEfficiencyModule.h:40
Belle2::KLMClusterEfficiencyModule::terminate
void terminate() override
This method is called at the end of the event processing.
Definition: KLMClusterEfficiencyModule.cc:339
Belle2::KLMClusterEfficiencyModule::endRun
void endRun() override
This method is called if the current run ends.
Definition: KLMClusterEfficiencyModule.cc:335
Belle2::StoreArray::getEntries
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:226
Belle2::BKLMHit2d
Store one BKLM strip hit as a ROOT object.
Definition: BKLMHit2d.h:42
Belle2::KLMClusterEfficiencyModule::m_DecayVertexZ
float m_DecayVertexZ
MCParticle decay vertex Z coordinate.
Definition: KLMClusterEfficiencyModule.h:109
Belle2::EKLMHit2d
Class for 2d hits handling.
Definition: EKLMHit2d.h:39
Belle2::KLMClusterEfficiencyModule::m_MCParticles
StoreArray< MCParticle > m_MCParticles
MC particles.
Definition: KLMClusterEfficiencyModule.h:164
Belle2::KLMClusterEfficiencyModule::event
void event() override
This method is called for each event.
Definition: KLMClusterEfficiencyModule.cc:97