Belle II Software development
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/KLMHit2d.h>
14
15/* Basf2 headers. */
16#include <framework/gearbox/Const.h>
17
18/* ROOT headers. */
19#include <Math/VectorUtil.h>
20#include <TCanvas.h>
21#include <TH2F.h>
22#include <TMarker.h>
23#include <TStyle.h>
24
25using namespace Belle2;
26
27REG_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");
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 ROOT::Math::XYZVector clusterPosition;
106 ROOT::Math::XYZVector decayVertex, hitPosition;
107 float angle;
108 /* cppcheck-suppress variableScope */
109 bool haveKL0;
110 n1 = m_KLMClusters.getEntries();
111 for (i1 = 0; i1 < n1; i1++) {
112 haveKL0 = false;
113 RelationVector<MCParticle> clusterMCParticles =
114 m_KLMClusters[i1]->getRelationsTo<MCParticle>();
115 n2 = clusterMCParticles.size();
116 for (i2 = 0; i2 < n2; i2++) {
117 if (clusterMCParticles[i2]->getPDG() == Const::Klong.getPDGCode())
118 haveKL0 = true;
119 }
120 if (haveKL0) {
121 if (n2 == 1)
123 else
125 } else
127 }
128 if (nevent < m_EventsClusterHistograms) {
129 static TH2F* hzx = new TH2F("hzx", "", 100, -300, 420, 100, -320, 320);
130 static TH2F* hzy = new TH2F("hzy", "", 100, -300, 420, 100, -320, 320);
131 static TH2F* hxy = new TH2F("hxy", "", 100, -320, 320, 100, -320, 320);
132 static TCanvas* c1 = new TCanvas();
133 static TMarker* clusterMarker = new TMarker(0, 0, 20);
134 static TMarker* hitMarker = new TMarker(0, 0, 21);
135 if (nevent == 0) {
136 hzx->GetXaxis()->SetTitle("z, cm");
137 hzx->GetYaxis()->SetTitle("x, cm");
138 hzy->GetXaxis()->SetTitle("z, cm");
139 hzy->GetYaxis()->SetTitle("y, cm");
140 hxy->GetXaxis()->SetTitle("x, cm");
141 hxy->GetYaxis()->SetTitle("y, cm");
142 }
143 hzx->Draw();
144 for (i1 = 0; i1 < n1; i1++) {
145 clusterMarker->SetMarkerColor(i1 + 1);
146 hitMarker->SetMarkerColor(i1 + 1);
147 clusterPosition = m_KLMClusters[i1]->getClusterPosition();
148 clusterMarker->DrawMarker(clusterPosition.Z(), clusterPosition.X());
149 RelationVector<KLMHit2d> klmHit2ds =
150 m_KLMClusters[i1]->getRelationsTo<KLMHit2d>();
151 n2 = klmHit2ds.size();
152 for (i2 = 0; i2 < n2; i2++) {
153 hitPosition = klmHit2ds[i2]->getPosition();
154 hitMarker->DrawMarker(hitPosition.Z(), hitPosition.X());
155 }
156 }
157 snprintf(str, 128, "clusters%dzx.eps", nevent);
158 c1->Print(str);
159 hzy->Draw();
160 for (i1 = 0; i1 < n1; i1++) {
161 clusterMarker->SetMarkerColor(i1 + 1);
162 hitMarker->SetMarkerColor(i1 + 1);
163 clusterPosition = m_KLMClusters[i1]->getClusterPosition();
164 clusterMarker->DrawMarker(clusterPosition.Z(), clusterPosition.Y());
165 RelationVector<KLMHit2d> klmHit2ds =
166 m_KLMClusters[i1]->getRelationsTo<KLMHit2d>();
167 n2 = klmHit2ds.size();
168 for (i2 = 0; i2 < n2; i2++) {
169 hitPosition = klmHit2ds[i2]->getPosition();
170 hitMarker->DrawMarker(hitPosition.Z(), hitPosition.Y());
171 }
172 }
173 snprintf(str, 128, "clusters%dzy.eps", nevent);
174 c1->Print(str);
175 hxy->Draw();
176 for (i1 = 0; i1 < n1; i1++) {
177 clusterMarker->SetMarkerColor(i1 + 1);
178 hitMarker->SetMarkerColor(i1 + 1);
179 clusterPosition = m_KLMClusters[i1]->getClusterPosition();
180 clusterMarker->DrawMarker(clusterPosition.X(), clusterPosition.Y());
181 RelationVector<KLMHit2d> klmHit2ds =
182 m_KLMClusters[i1]->getRelationsTo<KLMHit2d>();
183 n2 = klmHit2ds.size();
184 for (i2 = 0; i2 < n2; i2++) {
185 hitPosition = klmHit2ds[i2]->getPosition();
186 hitMarker->DrawMarker(hitPosition.X(), hitPosition.Y());
187 }
188 }
189 snprintf(str, 128, "clusters%dxy.eps", nevent);
190 c1->Print(str);
191 nevent++;
192 }
194 for (i1 = 0; i1 < n1; i1++) {
195 if (m_MCParticles[i1]->getPDG() != Const::Klong.getPDGCode())
196 continue;
197 decayVertex = m_MCParticles[i1]->getDecayVertex();
198 m_DecayVertexX = decayVertex.X();
199 m_DecayVertexY = decayVertex.Y();
200 m_DecayVertexZ = decayVertex.Z();
202 RelationVector<KLMHit2d> mcKLMHit2ds =
203 m_MCParticles[i1]->getRelationsFrom<KLMHit2d>();
204 n2 = mcKLMHit2ds.size();
205 for (i2 = 0; i2 < n2; i2++) {
206 hitPosition = mcKLMHit2ds[i2]->getPosition();
207 angle = ROOT::Math::VectorUtil::Angle(decayVertex, hitPosition);
208 if (angle > m_MaxDecayVertexHitAngle)
210 }
211 RelationVector<KLMCluster> kl0Clusters =
212 m_MCParticles[i1]->getRelationsFrom<KLMCluster>();
213 n2 = kl0Clusters.size();
214 if (n2 == 0)
216 else if (n2 == 1) {
217 RelationVector<KLMHit2d> klmHit2ds =
218 kl0Clusters[0]->getRelationsTo<KLMHit2d>();
219 bs1 = 0;
220 es1 = 0;
221 for (const KLMHit2d& hit2d : klmHit2ds) {
222 if (hit2d.getSubdetector() == KLMElementNumbers::c_BKLM)
223 bs1++;
224 else
225 es1++;
226 }
227 if (bs1 > 0) {
228 if (es1 > 0)
230 else
232 } else if (es1 > 0)
234 RelationVector<MCParticle> m_MCParticles2 =
235 kl0Clusters[0]->getRelationsTo<MCParticle>();
236 if (m_MCParticles2.size() == 1) {
237 if (m_MCParticles2.weight(0) == 1)
239 }
240 } else if (n2 == 2) {
241 RelationVector<KLMHit2d> klmHit2ds1 =
242 kl0Clusters[0]->getRelationsTo<KLMHit2d>();
243 RelationVector<KLMHit2d> klmHit2ds2 =
244 kl0Clusters[1]->getRelationsTo<KLMHit2d>();
245 bs1 = 0;
246 bs2 = 0;
247 es1 = 0;
248 es2 = 0;
249 for (const KLMHit2d& hit2d : klmHit2ds1) {
250 if (hit2d.getSubdetector() == KLMElementNumbers::c_BKLM)
251 bs1++;
252 else
253 es1++;
254 }
255 for (const KLMHit2d& hit2d : klmHit2ds2) {
256 if (hit2d.getSubdetector() == KLMElementNumbers::c_BKLM)
257 bs2++;
258 else
259 es2++;
260 }
261 if (bs1 > 0 && bs2 > 0) {
262 if (es1 > 0 && es2 > 0) {
264 } else if ((es1 > 0 && es2 == 0) || (es1 == 0 && es2 > 0)) {
266 } else if (es1 == 0 && es2 == 0) {
268 }
269 } else if (bs1 > 0 && bs2 == 0) {
270 if (es1 > 0 && es2 > 0) {
272 } else if (es1 == 0 && es2 > 0) {
274 }
275 } else if (bs1 == 0 && bs2 > 0) {
276 if (es1 > 0 && es2 > 0) {
278 } else if (es1 > 0 && es2 == 0) {
280 }
281 } else if (bs1 == 0 && bs2 == 0) {
282 if (es1 > 0 && es2 > 0) {
284 }
285 }
286 } else if (n2 >= 3)
288 if (m_SaveClusterData) {
289 for (i2 = 0; i2 < n2; i2++) {
290 clusterPosition = kl0Clusters[i2]->getClusterPosition();
291 m_ClusterX = clusterPosition.X();
292 m_ClusterY = clusterPosition.Y();
293 m_ClusterZ = clusterPosition.Z();
295 RelationVector<KLMHit2d> klmHit2ds =
296 kl0Clusters[i2]->getRelationsTo<KLMHit2d>();
297 n3 = klmHit2ds.size();
298 for (i3 = 0; i3 < n3; i3++) {
299 hitPosition = klmHit2ds[i3]->getPosition();
300 angle = ROOT::Math::VectorUtil::Angle(clusterPosition, hitPosition);
301 if (angle > m_MaxClusterHitAngle)
302 m_MaxClusterHitAngle = angle;
303 }
304 m_ClusterTree->Fill();
305 }
306 }
307 }
308}
309
311{
312}
313
315{
316 int i, rec1Cluster, rec2Clusters;
317 rec1Cluster = 0;
318 rec2Clusters = 0;
319 for (i = 0; i < 3; i++)
320 rec1Cluster = rec1Cluster + m_ReconstructedKL01Cluster[i];
321 for (i = 0; i < 6; i++)
322 rec2Clusters = rec2Clusters + m_ReconstructedKL02Clusters[i];
323 /* Always printed once, not necessary to use LogVar. */
324 B2INFO("Total number of KLM clusters: " << m_KL0Clusters +
326 B2INFO("K_L0 clusters: " << m_KL0Clusters);
327 B2INFO("(K_L0+other) clusters: " << m_PartlyKL0Clusters);
328 B2INFO("Other clusters: " << m_OtherClusters);
329 B2INFO("Total number of generated K_L0: " << m_NonreconstructedKL0 +
330 rec1Cluster + rec2Clusters + m_ReconstructedKL03Clusters);
331 B2INFO("Nonreconstructed K_L0: " << m_NonreconstructedKL0);
332 B2INFO("K_L0 reconstructed as 1 cluster (total): " << rec1Cluster);
333 B2INFO("K_L0 reconstructed as 1 cluster (BKLM): " <<
335 B2INFO("K_L0 reconstructed as 1 cluster (BKLM/EKLM): " <<
337 B2INFO("K_L0 reconstructed as 1 cluster (EKLM): " <<
339 B2INFO("K_L0 reconstructed as 1 cluster (exact reconstruction): " <<
341 B2INFO("K_L0 reconstructed as 2 clusters (total): " << rec2Clusters);
342 B2INFO("K_L0 reconstructed as 2 clusters (2 * BKLM): " <<
344 B2INFO("K_L0 reconstructed as 2 clusters (BKLM + BKLM/EKLM): " <<
346 B2INFO("K_L0 reconstructed as 2 clusters (BKLM + EKLM): " <<
348 B2INFO("K_L0 reconstructed as 2 clusters (2 * BKLM/EKLM): " <<
350 B2INFO("K_L0 reconstructed as 2 clusters (BKLM/EKLM + EKLM): " <<
352 B2INFO("K_L0 reconstructed as 2 clusters (2 * EKLM): " <<
354 B2INFO("K_L0 reconstructed as 3 or more clusters: " <<
357 m_ReconstructionTree->Fill();
358 m_OutputFile->cd();
359 m_ReconstructionTree->Write();
361 }
362 if (m_SaveClusterData) {
363 m_OutputFile->cd();
364 m_ClusterTree->Write();
365 delete m_ClusterTree;
366 }
368 delete m_OutputFile;
369}
370
int getPDGCode() const
PDG code.
Definition: Const.h:473
static const ParticleType Klong
K^0_L particle.
Definition: Const.h:678
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
KLM 2d hit.
Definition: KLMHit2d.h:33
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
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
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.