Belle II Software development
BKLMSimHistogrammerModule.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/bklm/modules/bklmSimHistogrammer/BKLMSimHistogrammerModule.h>
11
12/* Basf2 headers. */
13#include <framework/dataobjects/EventMetaData.h>
14#include <framework/datastore/StoreArray.h>
15#include <framework/datastore/StoreObjPtr.h>
16
17/* ROOT headers. */
18#include <TMath.h>
19
20/* C++ headers. */
21#include <iostream>
22
23using namespace Belle2;
24
25//-----------------------------------------------------------------
26// Register the Module
27//-----------------------------------------------------------------
28REG_MODULE(BKLMSimHistogrammer);
29
30//-----------------------------------------------------------------
31// Implementation
32//-----------------------------------------------------------------
33
35 m_hSimHitPerChannelLayer(nullptr),
36 m_hEvt(nullptr),
37 m_hSimHitPhiRPC(nullptr),
38 m_bgSourcePerLayer(nullptr),
39 m_bgSourcePerLayer2D(nullptr),
40 m_bgSourceVsPhi(nullptr),
41 m_bgSourceVsTheta(nullptr),
42 m_bgSource(nullptr),
43 m_bgSource2D(nullptr),
44 m_hSimHitPhiScinti(nullptr),
45 m_hSimHitThetaRPC(nullptr),
46 m_hSimHitThetaScinti(nullptr),
47 m_hSimHit_layer(nullptr),
48 m_hSimHit_layer2D(nullptr),
49 m_hSimHitThetaPhiRPC(nullptr),
50 m_hSimHitThetaPhiScinti(nullptr),
51 m_file(nullptr),
52 m_weight(0.0)
53{
54
55 setDescription("Analyzes bg");
56
57 addParam("filename", m_filename, "Output root filename", std::string("bg_output.root"));
58 addParam("realTime", m_realTime, "Time the simulation corresponds to", float(1.0));
59
60
61}
62
64{
65}
66
68{
69
70 hits2D.isRequired();
71 simHits.isRequired();
72
73
74 m_file = new TFile(m_filename.c_str(), "recreate");
75
76 m_hEvt = new TH1D("hEvts", "hEvts", 10001, -1, 10000);
77
78 m_hSimHitPhiRPC = new TH1D("simhitPhiRPC", "simHitPhiRPC", 100, -1 * TMath::Pi(), 1 * TMath::Pi());
79 m_hSimHitPhiScinti = new TH1D("simhitPhiScinti", "simHitPhiScinti", 100, -1 * TMath::Pi(), 1 * TMath::Pi());
80 m_hSimHitThetaRPC = new TH1D("simhitThetaRPC", "simHitThetaRPC", 100, 0, 1 * TMath::Pi());
81 m_hSimHitThetaScinti = new TH1D("simhitThetaScinti", "simHitThetaScinti", 100, 0, 1 * TMath::Pi());
82
83 m_hSimHit_layer = new TH1I("simHitLayer", "simHitLayer", 16, 0, 15);
84 m_hSimHit_layer2D = new TH1I("simHitLayer2D", "simHitLayer2D", 16, 0, 15);
85
86 m_bgSourcePerLayer = new TH2D("bgSourcPerLayer", "bgSourcePerLayer", 21, 0, 20, 16, 0, 15);
87 m_bgSourcePerLayer2D = new TH2D("bgSourcPerLayer2D", "bgSourcePerLayer2D", 21, 0, 20, 16, 0, 15);
88 m_bgSourceVsPhi = new TH2D("bgSourceVsPhi", "bgSourceVsPhi", 100, -1 * TMath::Pi(), TMath::Pi(), 21, 0, 20);
89 m_bgSourceVsTheta = new TH2D("bgSourceVsTheta", "bgSourceVsTheta", 100, 0, TMath::Pi(), 21, 0, 20);
90
91 m_bgSource = new TH1D("bgSource", "bgSource", 21, 0, 20);
92 m_bgSource2D = new TH1D("bgSource2D", "bgSource2D", 21, 0, 20);
93
94 m_hSimHitPerChannelLayer = new TH2D("posPerChannelLayer", "posPerChannelLayer", 5000, 0, 5000, 16, 0, 15);
95
96
97 // Force creation and persistence of BKLM output datastores
98}
99
101{
102 StoreObjPtr<EventMetaData> evtMetaData;
103 B2INFO("BKLMSimHistogrammer: Experiment " << evtMetaData->getExperiment() << " run " << evtMetaData->getRun());
104}
105
107{
108 //---------------------------------------------
109 // Get BKLM hits collection from the data store
110 //---------------------------------------------
111
112 int nSimHit = simHits.getEntries();
113
114 //use 1D hits (otherwise too many simHits, also look up background info), only problem might be that the threshold in simulation is not correct
115
116
117 // unsigned int d = 0;
118
119 if (m_realTime > 0)
120 m_weight = 1.0 / m_realTime;
121
122 std::cout << "real time " << m_realTime << " weight: " << m_weight << std::endl;
123
124 int n2DHits = hits2D.getEntries();
125 int n1DHits = hits1D.getEntries();
126 std::cout << "we have " << nSimHit << " sim hits " << n1DHits << " 1D hits " << n2DHits << " 2D hits " << std::endl;
127
128 // std::cout <<"we have " << digits.getEntries() <<" digits " <<std::endl;
129
130 m_hEvt->Fill(n2DHits, m_weight);
131 for (int i = 0; i < hits1D.getEntries(); i++) {
132 std::cout << "looking at 1DHit " << i << std::endl;
133 int scaledTag = -1;
134 RelationVector<KLMDigit> bklmDigits = hits1D[i]->getRelationsTo<KLMDigit>();
135 for (const auto& bklmDigit : bklmDigits) {
136 RelationVector<KLMSimHit> relatedSimHits = bklmDigit.getRelationsWith<KLMSimHit>();
137 for (const auto& simHit : relatedSimHits) {
138 auto bgTag = simHit.getBackgroundTag();
139 scaledTag = bgTag;
140 //other has numeric value of 99
141 if (scaledTag > 20)
142 scaledTag = 20;
143 std::cout << "scaledTag: " << scaledTag << std::endl;
144 break;
145 }
146 break;
147
148 }
149 int channel = hits1D[i]->getSection() * 840 + hits1D[i]->getSector() * 105 + hits1D[i]->isPhiReadout() * 1680 +
150 hits1D[i]->getStripAve();
151 m_hSimHitPerChannelLayer->Fill(channel, hits1D[i]->getLayer(), m_weight);
152 m_bgSource->Fill(scaledTag, m_weight);
153 m_bgSourcePerLayer->Fill(scaledTag, hits1D[i]->getLayer(), m_weight);
154 std::cout << "filling layer with tag: " << scaledTag << " and layer " << hits1D[i]->getLayer() << std::endl;
155 }
156
157 for (int i = 0; i < hits2D.getEntries(); i++) {
158 if (hits2D[i]->getSubdetector() != KLMElementNumbers::c_BKLM)
159 continue;
160 int scaledTag = -1;
161 ROOT::Math::XYZVector gHitPos = hits2D[i]->getPosition();
162 RelationVector<BKLMHit1d> related1DHits = hits2D[i]->getRelationsTo<BKLMHit1d>();
163 for (const auto& hit1d : related1DHits) {
164 RelationVector<KLMDigit> bklmDigits = hit1d.getRelationsTo<KLMDigit>();
165 for (const auto& bklmDigit : bklmDigits) {
166 RelationVector<KLMSimHit> relatedSimHits = bklmDigit.getRelationsWith<KLMSimHit>();
167 for (const auto& simHit : relatedSimHits) {
168 auto bgTag = simHit.getBackgroundTag();
169 scaledTag = bgTag;
170 //other has numeric value of 99
171 if (scaledTag > 20)
172 scaledTag = 20;
173
174
175 m_bgSourcePerLayer2D->Fill(scaledTag, hits2D[i]->getLayer(), m_weight);
176 m_bgSource2D->Fill(scaledTag, m_weight);
177 m_bgSourceVsTheta->Fill(gHitPos.Theta(), scaledTag, m_weight);
178 m_bgSourceVsPhi->Fill(gHitPos.Phi(), scaledTag, m_weight);
179 //one bg source is enough
180 break;
181 }
182 break;
183 }
184 break;
185 }
186 }
187 if (nSimHit == 0)
188 return;
189 for (int i = 0; i < n2DHits; i++) {
190 KLMHit2d* hit2D = hits2D[i];
191 ROOT::Math::XYZVector gHitPos = hit2D->getPosition();
192 if (hit2D->inRPC()) {
193 m_hSimHitPhiRPC->Fill(gHitPos.Phi(), m_weight);
194 m_hSimHitThetaRPC->Fill(gHitPos.Theta(), m_weight);
195 } else {
196 m_hSimHitPhiScinti->Fill(gHitPos.Phi(), m_weight);
197 m_hSimHitThetaScinti->Fill(gHitPos.Theta(), m_weight);
198 }
199 }
200 for (int h = 0; h < nSimHit; ++h) {
201 KLMSimHit* simHit = simHits[h];
202
203 RelationVector<MCParticle> bklmMCParticles = simHit->getRelationsFrom<MCParticle>();
204 if (bklmMCParticles.size() > 0) {
205 std::cout << "found MC particle!" << std::endl;
206 RelationVector<MCParticle> bklmMCParticlesTo = simHit->getRelationsTo<MCParticle>();
207 std::cout << "got " << bklmMCParticlesTo.size() << " as relation to " << std::endl;
208 }
209 }
210}
211
213{
214}
215
217{
218 m_file->Write();
219 m_file->Close();
220}
Store one reconstructed BKLM 1D hit as a ROOT object.
Definition: BKLMHit1d.h:30
TH1D * m_hSimHitPhiScinti
histogram for sim hit phi
float m_realTime
time this simulation corresponds to
TH2D * m_bgSourcePerLayer2D
bg source of 2D hits
TH1I * m_hSimHit_layer2D
histogram for the layers of 2D hits
void initialize() override
Initialize at start of job.
void event() override
Digitize one event and write hits, digis, and relations into DataStore.
TH1D * m_hSimHitThetaRPC
histogram for sim hit theta
void endRun() override
Do any needed actions at the end of a simulation run.
void terminate() override
Terminate at the end of job.
TH1D * m_hSimHitThetaScinti
histogram for sim hit theta scint
StoreArray< KLMSimHit > simHits
KLMSimHit StoreArray.
TH2D * m_bgSourceVsPhi
bg source of 2D hits vs phi
StoreArray< BKLMHit1d > hits1D
hits1D StoreArray
TH1I * m_hSimHit_layer
histogram for the layers of 1D hits
void beginRun() override
Do any needed actions at the start of a simulation run.
TH1D * m_hSimHitPhiRPC
histogram for sim hit phi
StoreArray< KLMHit2d > hits2D
hits2D StoreArray
TFile * m_file
Digitize hit(s) in one scintillator strip with pulse-shape fit.
TH2D * m_bgSourceVsTheta
bg source of 2D hits vs theta
float m_weight
weight for each event (inverse of the realTime)
std::string m_filename
filename for the root file
TH2D * m_bgSourcePerLayer
bg source of 1D hits
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition: KLMDigit.h:29
KLM 2d hit.
Definition: KLMHit2d.h:33
bool inRPC() const
Determine whether this 2D hit is in RPC or scintillator.
Definition: KLMHit2d.h:185
ROOT::Math::XYZVector getPosition() const
Get hit global position.
Definition: KLMHit2d.h:315
KLM simulation hit.
Definition: KLMSimHit.h:31
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.
RelationVector< FROM > getRelationsFrom(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from another store array to this object.
RelationVector< TO > getRelationsTo(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from this object to another store array.
virtual unsigned short getBackgroundTag() const
Get background tag.
Definition: SimHitBase.h:46
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
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.