Belle II Software  release-06-02-00
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 /* KLM headers. */
13 #include <klm/dataobjects/bklm/BKLMSimHit.h>
14 #include <klm/dataobjects/KLMDigit.h>
15 
16 /* Belle 2 headers. */
17 #include <framework/dataobjects/EventMetaData.h>
18 #include <framework/datastore/StoreArray.h>
19 #include <framework/datastore/StoreObjPtr.h>
20 
21 /* ROOT headers. */
22 #include <TMath.h>
23 
24 /* C++ headers. */
25 #include <iostream>
26 
27 using namespace std;
28 using namespace Belle2;
29 
30 //-----------------------------------------------------------------
31 // Register the Module
32 //-----------------------------------------------------------------
33 REG_MODULE(BKLMSimHistogrammer)
34 
35 //-----------------------------------------------------------------
36 // Implementation
37 //-----------------------------------------------------------------
38 
40  m_hSimHitPerChannelLayer(nullptr),
41  m_hEvt(nullptr),
42  m_hSimHitPhiRPC(nullptr),
43  m_bgSourcePerLayer(nullptr),
44  m_bgSourcePerLayer2D(nullptr),
45  m_bgSourceVsPhi(nullptr),
46  m_bgSourceVsTheta(nullptr),
47  m_bgSource(nullptr),
48  m_bgSource2D(nullptr),
49  m_hSimHitPhiScinti(nullptr),
50  m_hSimHitThetaRPC(nullptr),
51  m_hSimHitThetaScinti(nullptr),
52  m_hSimHit_layer(nullptr),
53  m_hSimHit_layer2D(nullptr),
54  m_hSimHitThetaPhiRPC(nullptr),
55  m_hSimHitThetaPhiScinti(nullptr),
56  m_file(nullptr),
57  m_weight(0.0)
58 {
59 
60  setDescription("Analyzes bg");
61 
62  addParam("filename", m_filename, "Output root filename", string("bg_output.root"));
63  addParam("realTime", m_realTime, "Time the simulation corresponds to", float(1.0));
64 
65 
66 }
67 
68 BKLMSimHistogrammerModule::~BKLMSimHistogrammerModule()
69 {
70 }
71 
72 void BKLMSimHistogrammerModule::initialize()
73 {
74 
75  hits2D.isRequired();
76  simHits.isRequired();
77 
78 
79  m_file = new TFile(m_filename.c_str(), "recreate");
80 
81  m_hEvt = new TH1D("hEvts", "hEvts", 10001, -1, 10000);
82 
83  m_hSimHitPhiRPC = new TH1D("simhitPhiRPC", "simHitPhiRPC", 100, -1 * TMath::Pi(), 1 * TMath::Pi());
84  m_hSimHitPhiScinti = new TH1D("simhitPhiScinti", "simHitPhiScinti", 100, -1 * TMath::Pi(), 1 * TMath::Pi());
85  m_hSimHitThetaRPC = new TH1D("simhitThetaRPC", "simHitThetaRPC", 100, 0, 1 * TMath::Pi());
86  m_hSimHitThetaScinti = new TH1D("simhitThetaScinti", "simHitThetaScinti", 100, 0, 1 * TMath::Pi());
87 
88  m_hSimHit_layer = new TH1I("simHitLayer", "simHitLayer", 16, 0, 15);
89  m_hSimHit_layer2D = new TH1I("simHitLayer2D", "simHitLayer2D", 16, 0, 15);
90 
91  m_bgSourcePerLayer = new TH2D("bgSourcPerLayer", "bgSourcePerLayer", 21, 0, 20, 16, 0, 15);
92  m_bgSourcePerLayer2D = new TH2D("bgSourcPerLayer2D", "bgSourcePerLayer2D", 21, 0, 20, 16, 0, 15);
93  m_bgSourceVsPhi = new TH2D("bgSourceVsPhi", "bgSourceVsPhi", 100, -1 * TMath::Pi(), TMath::Pi(), 21, 0, 20);
94  m_bgSourceVsTheta = new TH2D("bgSourceVsTheta", "bgSourceVsTheta", 100, 0, TMath::Pi(), 21, 0, 20);
95 
96  m_bgSource = new TH1D("bgSource", "bgSource", 21, 0, 20);
97  m_bgSource2D = new TH1D("bgSource2D", "bgSource2D", 21, 0, 20);
98 
99  m_hSimHitPerChannelLayer = new TH2D("posPerChannelLayer", "posPerChannelLayer", 5000, 0, 5000, 16, 0, 15);
100 
101 
102  // Force creation and persistence of BKLM output datastores
103 }
104 
105 void BKLMSimHistogrammerModule::beginRun()
106 {
107  StoreObjPtr<EventMetaData> evtMetaData;
108  B2INFO("BKLMSimHistogrammer: Experiment " << evtMetaData->getExperiment() << " run " << evtMetaData->getRun());
109 }
110 
111 void BKLMSimHistogrammerModule::event()
112 {
113  //---------------------------------------------
114  // Get BKLM hits collection from the data store
115  //---------------------------------------------
116 
117  int nSimHit = simHits.getEntries();
118 
119  //use 1D hits (otherwise too many simHits, also look up background info), only problem might be that the threshold in simulation is not correct
120 
121 
122  // unsigned int d = 0;
123 
124  if (m_realTime > 0)
125  m_weight = 1.0 / m_realTime;
126 
127  cout << "real time " << m_realTime << " weight: " << m_weight << endl;
128 
129  int n2DHits = hits2D.getEntries();
130  int n1DHits = hits1D.getEntries();
131  cout << "we have " << nSimHit << " sim hits " << n1DHits << " 1D hits " << n2DHits << " 2D hits " << endl;
132 
133  // cout <<"we have " << digits.getEntries() <<" digits " <<endl;
134 
135  m_hEvt->Fill(n2DHits, m_weight);
136  for (int i = 0; i < hits1D.getEntries(); i++) {
137  cout << "looking at 1DHit " << i << endl;
138  int scaledTag = -1;
139  RelationVector<KLMDigit> bklmDigits = hits1D[i]->getRelationsTo<KLMDigit>();
140  for (const auto& bklmDigit : bklmDigits) {
141  RelationVector<BKLMSimHit> relatedSimHits = bklmDigit.getRelationsWith<BKLMSimHit>();
142  for (const auto& simHit : relatedSimHits) {
143  auto bgTag = simHit.getBackgroundTag();
144  scaledTag = bgTag;
145  //other has numeric value of 99
146  if (scaledTag > 20)
147  scaledTag = 20;
148  cout << "scaledTag: " << scaledTag << endl;
149  break;
150 
151 
152 
153 // switch(bgTag)
154 // {
155 // case SimHitBase::bg_none:
156 // cout <<"this is no bg " <<endl;
157 // break;
158 // case SimHitBase::bg_Coulomb_LER:
159 // cout <<"coulomb ler" <<endl;
160 // break;
161 // case SimHitBase::bg_Coulomb_HER:
162 // cout <<"coulomb her" <<endl;
163 // break;
164 //
165 // case SimHitBase::bg_Touschek_LER:
166 // cout <<"touschek ler" <<endl;
167 // break;
168 // case SimHitBase::bg_Touschek_HER:
169 // cout <<"touschek her" <<endl;
170 // break;
171 //
172 // case SimHitBase::bg_RBB_LER:
173 // cout <<"rbb ler" <<endl;
174 // break;
175 //
176 // case SimHitBase::bg_RBB_HER:
177 // cout <<"rbb her" <<endl;
178 // break;
179 //
180 // case SimHitBase::bg_twoPhoton:
181 // cout <<" two phot" <<endl;
182 // break;
183 // case SimHitBase::bg_RBB_gamma:
184 // cout <<"rbb gamma" <<endl;
185 // break;
186 // case SimHitBase::bg_RBB_LER_far:
187 // cout <<"rbb ler far" <<endl;
188 // break;
189 // case SimHitBase::bg_RBB_HER_far:
190 // cout <<"rbb her far" <<endl;
191 // break;
192 // case SimHitBase::bg_Touschek_LER_far:
193 // cout <<"touschek ler far" <<endl;
194 // break;
195 // case SimHitBase::bg_Touschek_HER_far:
196 // cout <<"touschek her far" <<endl;
197 // break;
198 //
199 // case SimHitBase::bg_SynchRad_LER:
200 // cout <<"synchRad ler " <<endl;
201 // break;
202 // case SimHitBase::bg_SynchRad_HER:
203 // cout <<"synchRad her " <<endl;
204 // break;
205 // case SimHitBase::bg_BHWide_LER:
206 // cout <<"bh wide ler " <<endl;
207 // break;
208 // case SimHitBase::bg_BHWide_HER:
209 // cout <<"bh wide her " <<endl;
210 // break;
211 // case SimHitBase::bg_other:
212 // cout <<"bg other " <<endl;
213 // break;
214 //
215 // default:
216 // {
217 // cout <<"this bg code is not allowed! " << bgTag <<endl;
218 // break;
219 // }
220 //
221 // }
222  }
223  break;
224 
225  }
226  int channel = hits1D[i]->getSection() * 840 + hits1D[i]->getSector() * 105 + hits1D[i]->isPhiReadout() * 1680 +
227  hits1D[i]->getStripAve();
228  m_hSimHitPerChannelLayer->Fill(channel, hits1D[i]->getLayer(), m_weight);
229  m_bgSource->Fill(scaledTag, m_weight);
230  m_bgSourcePerLayer->Fill(scaledTag, hits1D[i]->getLayer(), m_weight);
231  cout << "filling layer with tag: " << scaledTag << " and layer " << hits1D[i]->getLayer() << endl;
232  }
233 
234  for (int i = 0; i < hits2D.getEntries(); i++) {
235  int scaledTag = -1;
236  TVector3 gHitPos = hits2D[i]->getGlobalPosition();
237  RelationVector<BKLMHit1d> related1DHits = hits2D[i]->getRelationsTo<BKLMHit1d>();
238  for (const auto& hit1d : related1DHits) {
239  RelationVector<KLMDigit> bklmDigits = hit1d.getRelationsTo<KLMDigit>();
240  for (const auto& bklmDigit : bklmDigits) {
241  RelationVector<BKLMSimHit> relatedSimHits = bklmDigit.getRelationsWith<BKLMSimHit>();
242  for (const auto& simHit : relatedSimHits) {
243  auto bgTag = simHit.getBackgroundTag();
244  scaledTag = bgTag;
245  //other has numeric value of 99
246  if (scaledTag > 20)
247  scaledTag = 20;
248 
249 
250  m_bgSourcePerLayer2D->Fill(scaledTag, hits2D[i]->getLayer(), m_weight);
251  m_bgSource2D->Fill(scaledTag, m_weight);
252  m_bgSourceVsTheta->Fill(gHitPos.Theta(), scaledTag, m_weight);
253  m_bgSourceVsPhi->Fill(gHitPos.Phi(), scaledTag, m_weight);
254  //one bg source is enough
255  break;
256  }
257  break;
258  }
259  break;
260  }
261  }
262 
263 
264 
265 
266 // for(int i=0;i<digits.getEntries();i++)
267 // {
268 // RelationVector<BKLMSimHit> simHits=digits[i]->getRelationsFrom<BKLMSimHit>();
269 // cout <<"digit " << i <<" has " << simHits.size() << " sim Hits associated " <<endl;
270 // }
271 // // RelationVector<BKLMSimHit> bklmSimHits = bklmDigits[i4]->getRelationsFrom<BKLMSimHit>();
272 //
273 // for (int i3 = 0; i3 < n1DHits; i3++) {
274 // RelationVector<BKLMDigit> bklmDigits =
275 // hits1D[i3]->getRelationsTo<BKLMDigit>();
276 // cout <<"the 1Dhit " << i3 << "has "<< bklmDigits.size()<< " digits" <<endl;
277 // int n4 = bklmDigits.size();
278 // for (int i4 = 0; i4 < n4; i4++) {
279 // RelationVector<BKLMSimHit> bklmSimHits =
280 // bklmDigits[i4]->getRelationsFrom<BKLMSimHit>();
281 // cout <<"and this digit is associated with " << bklmSimHits.size()<<" sim hits " <<endl;
282 // // n5 = bklmSimHits.size();
283 // }}
284 //
285 //
286 
287 
288  if (nSimHit == 0)
289  return;
290  for (int i = 0; i < n2DHits; i++) {
291  BKLMHit2d* hit2D = hits2D[i];
292  TVector3 gHitPos = hit2D->getGlobalPosition();
293  if (hit2D->inRPC()) {
294  m_hSimHitPhiRPC->Fill(gHitPos.Phi(), m_weight);
295  m_hSimHitThetaRPC->Fill(gHitPos.Theta(), m_weight);
296  } else {
297  m_hSimHitPhiScinti->Fill(gHitPos.Phi(), m_weight);
298  m_hSimHitThetaScinti->Fill(gHitPos.Theta(), m_weight);
299  }
300 
301 
302 
303  }
304 
305  for (int h = 0; h < nSimHit; ++h) {
306  BKLMSimHit* simHit = simHits[h];
307 
308  RelationVector<MCParticle> bklmMCParticles = simHit->getRelationsFrom<MCParticle>();
309 // cout <<"did " << bklmMCParticles.size() <<" as relation from " <<endl;
310  if (bklmMCParticles.size() > 0) {
311  cout << "found MC particle!" << endl;
312  RelationVector<MCParticle> bklmMCParticlesTo = simHit->getRelationsTo<MCParticle>();
313  cout << "got " << bklmMCParticlesTo.size() << " as relation to " << endl;
314  }
315 
316 
317 
318 // //getMomentum, and getProductionVector getPDG
319  }
320 
321 
322 //
323 //td::map<int, std::vector<std::pair<int, BKLMSimHit*> > > volIDToSimHits;
324 //
325 //
326 // BKLMSimHit* simHit = simHits[h];
327 //
328 // RelationVector<MCParticle> bklmMCParticles =simHit->getRelationsFrom<MCParticle>();
329 // cout <<"did " << bklmMCParticles.size() <<" as relation from " <<endl;
330 //
331 // RelationVector<MCParticle> bklmMCParticlesTo =simHit->getRelationsTo<MCParticle>();
332 // // RelationVector<MCParticle> bklmMCParticlesTo =simHit->getRelationsTo<MCParticle>();
333 // cout <<"got " << bklmMCParticlesTo.size() <<" as relation to " <<endl;
334 // //getMomentum, and getProductionVector getPDG
335 //
336 //
337 // // cout <<" looking at sector: "<< simHit->getSector() << " get layer: " << simHit->getLayer()<<" strip: "<< simHit->getStrip()<<endl;
338 // m_hSimHit_layer->Fill(simHit->getLayer());
339 // // if (simHit->inRPC()) {
340 // {
341 // //get relations
342 // RelationVector<BKLMSimHitPosition> simPositionRelation=simHit->getRelationsTo<BKLMSimHitPosition>();
343 // cout <<" got " << simPositionRelation.size() <<" pos" <<endl;
344 //
345 //
346 //
347 // for(unsigned int iPos=0;iPos<simPositionRelation.size();iPos++)
348 // {
349 // BKLMSimHitPosition *simPos=simPositionRelation[iPos];
350 // if(simPos)
351 // cout <<"found simPos: "<< simPos->getPosition().Phi() <<" theta: "<< simPos->getPosition().Theta()<<endl;
352 // else
353 // cout <<"no sim pos " <<endl;
354 // }
355 //
356 // BKLMSimHitPosition* pos = simHit->getRelatedFrom<BKLMSimHitPosition>();
357 // if (!pos) {
358 // //nothing found, do some error handling here
359 // cout <<"no position found .." <<endl;
360 // }
361 // else
362 // {
363 // cout <<" found pos " << pos->getPosition().Phi() <<" theta: "<< pos->getPosition().Theta() <<endl;
364 // }
365 //
366 //
367 // }
368 //
369 
370 }
371 
372 void BKLMSimHistogrammerModule::endRun()
373 {
374 }
375 
376 void BKLMSimHistogrammerModule::terminate()
377 {
378  m_file->Write();
379  m_file->Close();
380 }
Store one reconstructed BKLM 1D hit as a ROOT object.
Definition: BKLMHit1d.h:30
Store one BKLM strip hit as a ROOT object.
Definition: BKLMHit2d.h:32
bool inRPC() const
Determine whether this 2D hit is in RPC or scintillator.
Definition: BKLMHit2d.h:59
TVector3 getGlobalPosition(void) const
Get 3D hit position in global coordinates.
Definition: BKLMHit2d.h:179
Convert BKLM raw simulation hits to digitizations.
Store one simulation hit as a ROOT object.
Definition: BKLMSimHit.h:25
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition: KLMDigit.h:30
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.
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:95
#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.