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