Belle II Software development
MCMatcherKLMClustersModule.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/MCMatcherKLMClusters/MCMatcherKLMClustersModule.h>
11
12/* KLM headers. */
13#include <klm/dataobjects/bklm/BKLMHit1d.h>
14#include <klm/dataobjects/KLMDigit.h>
15#include <klm/dataobjects/KLMHit2d.h>
16#include <klm/dataobjects/KLMSimHit.h>
17
18/* Basf2 headers. */
19#include <mdst/dataobjects/MCParticle.h>
20
21/* C++ headers. */
22#include <map>
23
24using namespace Belle2;
25
26REG_MODULE(MCMatcherKLMClusters);
27
29{
30 setDescription("Module for MC matching for KLM clusters.");
32 addParam("Hit2dRelations", m_Hit2dRelations,
33 "Add also relations for KLMHit2d and KLMHit2d.", false);
34}
35
37{
38}
39
41{
42 m_KLMClusters.isRequired();
44 if (m_MCParticles.isValid()) {
45 m_KLMClusters.registerRelationTo(m_MCParticles);
46 }
47 if (m_Hit2dRelations) {
48 StoreArray<KLMHit2d> bklmHit2ds;
49 StoreArray<KLMHit2d> eklmHit2ds;
50 if (m_MCParticles.isValid()) {
53 }
54 }
55}
56
58{
59 // Don't do anything if MCParticles aren't present
60 if (not m_MCParticles.isValid()) {
61 return;
62 }
63
64 double weightSum;
65 /* cppcheck-suppress variableScope */
66 int i1, i2, i3, i4, i5, i6, n1, n2, n3, n4, n5, n6;
67 std::map<MCParticle*, double> mcParticles, mcParticlesHit;
68 std::map<MCParticle*, double>::iterator it;
69 n1 = m_KLMClusters.getEntries();
70 for (i1 = 0; i1 < n1; i1++) {
71 mcParticles.clear();
72 RelationVector<KLMHit2d> bklmHit2ds =
73 m_KLMClusters[i1]->getRelationsTo<KLMHit2d>();
74 n2 = bklmHit2ds.size();
75 for (i2 = 0; i2 < n2; i2++) {
77 mcParticlesHit.clear();
78 RelationVector<BKLMHit1d> bklmHit1ds =
79 bklmHit2ds[i2]->getRelationsTo<BKLMHit1d>();
80 n3 = bklmHit1ds.size();
81 for (i3 = 0; i3 < n3; i3++) {
82 RelationVector<KLMDigit> bklmDigits =
83 bklmHit1ds[i3]->getRelationsTo<KLMDigit>();
84 n4 = bklmDigits.size();
85 for (i4 = 0; i4 < n4; i4++) {
86 RelationVector<KLMSimHit> bklmSimHits =
87 bklmDigits[i4]->getRelationsTo<KLMSimHit>();
88 n5 = bklmSimHits.size();
89 for (i5 = 0; i5 < n5; i5++) {
90 RelationVector<MCParticle> bklmMCParticles =
91 bklmSimHits[i5]->getRelationsFrom<MCParticle>();
92 n6 = bklmMCParticles.size();
93 for (i6 = 0; i6 < n6; i6++) {
94 it = mcParticles.find(bklmMCParticles[i6]);
95 if (it == mcParticles.end()) {
96 mcParticles.insert(std::pair<MCParticle*, double>(
97 bklmMCParticles[i6],
98 bklmSimHits[i5]->getEnergyDeposit()));
99 } else {
100 it->second = it->second + bklmSimHits[i5]->getEnergyDeposit();
101 }
102 if (m_Hit2dRelations) {
103 it = mcParticlesHit.find(bklmMCParticles[i6]);
104 if (it == mcParticlesHit.end()) {
105 mcParticlesHit.insert(std::pair<MCParticle*, double>(
106 bklmMCParticles[i6],
107 bklmSimHits[i5]->getEnergyDeposit()));
108 } else {
109 it->second = it->second + bklmSimHits[i5]->getEnergyDeposit();
110 }
111 }
112 }
113 }
114 }
115 }
116 if (m_Hit2dRelations) {
117 weightSum = 0;
118 for (it = mcParticlesHit.begin(); it != mcParticlesHit.end(); ++it)
119 weightSum = weightSum + it->second;
120 for (it = mcParticlesHit.begin(); it != mcParticlesHit.end(); ++it)
121 bklmHit2ds[i2]->addRelationTo(it->first, it->second / weightSum);
122 }
123 }
124 RelationVector<KLMHit2d> eklmHit2ds =
125 m_KLMClusters[i1]->getRelationsTo<KLMHit2d>();
126 n2 = eklmHit2ds.size();
127 for (i2 = 0; i2 < n2; i2++) {
129 mcParticlesHit.clear();
130 RelationVector<KLMDigit> eklmDigits =
131 eklmHit2ds[i2]->getRelationsTo<KLMDigit>();
132 n3 = eklmDigits.size();
133 for (i3 = 0; i3 < n3; i3++) {
134 RelationVector<KLMSimHit> eklmSimHits =
135 eklmDigits[i3]->getRelationsTo<KLMSimHit>();
136 n4 = eklmSimHits.size();
137 for (i4 = 0; i4 < n4; i4++) {
138 RelationVector<MCParticle> eklmMCParticles =
139 eklmSimHits[i4]->getRelationsFrom<MCParticle>();
140 n5 = eklmMCParticles.size();
141 for (i5 = 0; i5 < n5; i5++) {
142 it = mcParticles.find(eklmMCParticles[i5]);
143 if (it == mcParticles.end()) {
144 mcParticles.insert(std::pair<MCParticle*, double>(
145 eklmMCParticles[i5],
146 eklmSimHits[i4]->getEnergyDeposit()));
147 } else {
148 it->second = it->second + eklmSimHits[i4]->getEnergyDeposit();
149 }
150 if (m_Hit2dRelations) {
151 it = mcParticlesHit.find(eklmMCParticles[i5]);
152 if (it == mcParticlesHit.end()) {
153 mcParticlesHit.insert(std::pair<MCParticle*, double>(
154 eklmMCParticles[i5],
155 eklmSimHits[i4]->getEnergyDeposit()));
156 } else {
157 it->second = it->second + eklmSimHits[i4]->getEnergyDeposit();
158 }
159 }
160 }
161 }
162 }
163 if (m_Hit2dRelations) {
164 weightSum = 0;
165 for (it = mcParticlesHit.begin(); it != mcParticlesHit.end(); ++it)
166 weightSum = weightSum + it->second;
167 for (it = mcParticlesHit.begin(); it != mcParticlesHit.end(); ++it)
168 eklmHit2ds[i2]->addRelationTo(it->first, it->second / weightSum);
169 }
170 }
171 weightSum = 0;
172 for (it = mcParticles.begin(); it != mcParticles.end(); ++it)
173 weightSum = weightSum + it->second;
174 for (it = mcParticles.begin(); it != mcParticles.end(); ++it)
175 m_KLMClusters[i1]->addRelationTo(it->first, it->second / weightSum);
176 }
177}
Store one reconstructed BKLM 1D hit as a ROOT object.
Definition: BKLMHit1d.h:30
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition: KLMDigit.h:29
KLM 2d hit.
Definition: KLMHit2d.h:33
KLM simulation hit.
Definition: KLMSimHit.h:31
bool m_Hit2dRelations
Add relations for KLMHit2d and KLMHit2d.
void event() override
This method is called for each event.
StoreArray< KLMCluster > m_KLMClusters
KLM clusters.
StoreArray< MCParticle > m_MCParticles
MCParticles StoreArray.
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
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ 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:80
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
bool isValid() const
Check wether the array was registered.
Definition: StoreArray.h:288
bool registerRelationTo(const StoreArray< TO > &toArray, DataStore::EDurability durability=DataStore::c_Event, DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut, const std::string &namedRelation="") const
Register a relation to the given StoreArray.
Definition: StoreArray.h:140
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.