Belle II Software  release-05-02-19
ClawsStudyModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2013 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Igal Jaegle *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <beast/claws/modules/ClawsStudyModule.h>
12 #include <beast/claws/dataobjects/CLAWSSimHit.h>
13 #include <beast/claws/dataobjects/ClawsHit.h>
14 #include <generators/SAD/dataobjects/SADMetaHit.h>
15 #include <framework/datastore/StoreArray.h>
16 #include <framework/gearbox/GearDir.h>
17 #include <framework/logging/Logger.h>
18 
19 #include <fstream>
20 #include <string>
21 
22 // ROOT
23 #include <TH1.h>
24 #include <TH2.h>
25 
26 using namespace std;
27 
28 using namespace Belle2;
29 using namespace claws;
30 
31 //-----------------------------------------------------------------
32 // Register the Module
33 //-----------------------------------------------------------------
34 REG_MODULE(ClawsStudy)
35 
36 //-----------------------------------------------------------------
37 // Implementation
38 //-----------------------------------------------------------------
39 
41 {
42  // Set module properties
43  setDescription("Study module for Clawss (BEAST)");
44 
45  addParam("Ethres", m_Ethres, "Energy threshold in MeV", 0.0);
46 }
47 
48 ClawsStudyModule::~ClawsStudyModule()
49 {
50 }
51 
52 //This module is a histomodule. Any histogram created here will be saved by the HistoManager module
53 void ClawsStudyModule::defineHisto()
54 {
55  /*
56  for (int i = 0; i < 16; i++) {
57  h_clawss_Evtof1[i] = new TH2F(TString::Format("clawss_Evtof1_%d", i), "Energy deposited [MeV] vs TOF [ns] - all", 500, 0., 1000.,
58  100, 0., 10.);
59  h_clawss_Evtof2[i] = new TH2F(TString::Format("clawss_Evtof2_%d", i), "Energy deposited [MeV] vs TOF [ns] - only photons", 500, 0.,
60  100., 1000, 0., 10.);
61  h_clawss_Evtof3[i] = new TH2F(TString::Format("clawss_Evtof3_%d", i), "Energy deposited [MeV] vs TOF [ns] - only e+/e-", 500, 0.,
62  100., 1000, 0., 10.);
63  h_clawss_Evtof4[i] = new TH2F(TString::Format("clawss_Evtof4_%d", i), "Energy deposited [MeV] vs TOF [ns] - only e+/e-", 500, 0.,
64  100., 1000, 0., 10.);
65  h_clawss_edep[i] = new TH1F(TString::Format("clawss_edep_%d", i), "Energy deposited [MeV]", 5000, 0., 10.);
66  h_Wclawss_edep[i] = new TH1F(TString::Format("Wclawss_edep_%d", i), "Energy deposited [MeV]", 5000, 0., 10.);
67  }
68  */
69  h_clawss_hitrate1 = new TH1F("clawss_hitrate1", "Hit distributions", 16, 0., 16.);
70  h_clawss_hitrate2 = new TH1F("clawss_hitrate2", "Hit distributions", 16, 0., 16.);
71  h_clawss_hitrate1W = new TH1F("clawss_hitrate1W", "Hit distributions", 16, 0., 16.);
72  h_clawss_hitrate2W = new TH1F("clawss_hitrate2W", "Hit distributions", 16, 0., 16.);
73 
74  h_clawss_hitrate1->Sumw2();
75  h_clawss_hitrate1W->Sumw2();
76  h_clawss_hitrate2->Sumw2();
77  h_clawss_hitrate2W->Sumw2();
78 
79  h_clawss_rs_hitrate1 = new TH2F("clawss_rs_hitrate1", "Hit distributions vs rs", 16, 0., 16., 12, 0., 12.);
80  h_clawss_rs_hitrate2 = new TH2F("clawss_rs_hitrate2", "Hit distributions vs rs", 16, 0., 16., 12, 0., 12.);
81  h_clawss_rs_hitrate1W = new TH2F("clawss_rs_hitrate1W", "Hit distributions vs rs", 16, 0., 16., 12, 0., 12.);
82  h_clawss_rs_hitrate2W = new TH2F("clawss_rs_hitrate2W", "Hit distributions vs rs", 16, 0., 16., 12, 0., 12.);
83 
84  h_clawss_rs_hitrate1->Sumw2();
85  h_clawss_rs_hitrate1W->Sumw2();
86  h_clawss_rs_hitrate2->Sumw2();
87  h_clawss_rs_hitrate2W->Sumw2();
88 
89  for (int i = 0; i < 16; i++) {
90  h_clawss_rate1[i] = new TH1F(TString::Format("clawss_rate1_%d", i), "PE distributions", 500, 0., 500.);
91  h_clawss_rate2[i] = new TH1F(TString::Format("clawss_rate2_%d", i), "PE distributions", 500, 0., 500.);
92  h_clawss_rate1W[i] = new TH1F(TString::Format("clawss_rate1W_%d", i), "PE distributions", 500, 0., 500.);
93  h_clawss_rate2W[i] = new TH1F(TString::Format("clawss_rate2W_%d", i), "PE distributions", 500, 0., 500.);
94  h_clawss_pe1[i] = new TH2F(TString::Format("clawss_pe1_%d", i), "PE distributions", 500, 0., 500., 100, 0., 1000.);
95  h_clawss_pe2[i] = new TH2F(TString::Format("clawss_pe2_%d", i), "PE distributions", 500, 0., 500., 100, 0., 1000.);
96  h_clawss_pe1W[i] = new TH2F(TString::Format("clawss_pe1W_%d", i), "PE distributions", 500, 0., 500., 100, 0., 1000.);
97  h_clawss_pe2W[i] = new TH2F(TString::Format("clawss_pe2W_%d", i), "PE distributions", 500, 0., 500., 100, 0., 1000.);
98 
99  h_clawss_rs_rate1[i] = new TH2F(TString::Format("clawss_rs_rate1_%d", i), "PE distributions", 500, 0., 500., 12, 0., 12.);
100  h_clawss_rs_rate2[i] = new TH2F(TString::Format("clawss_rs_rate2_%d", i), "PE distributions", 500, 0., 500., 12, 0., 12.);
101  h_clawss_rs_rate1W[i] = new TH2F(TString::Format("clawss_rs_rate1W_%d", i), "PE distributions", 500, 0., 500., 12, 0., 12.);
102  h_clawss_rs_rate2W[i] = new TH2F(TString::Format("clawss_rs_rate2W_%d", i), "PE distributions", 500, 0., 500., 12, 0., 12.);
103 
104  h_clawss_rate1[i]->Sumw2();
105  h_clawss_rate2[i]->Sumw2();
106  h_clawss_rate1W[i]->Sumw2();
107  h_clawss_rate2W[i]->Sumw2();
108  h_clawss_rs_rate1[i]->Sumw2();
109  h_clawss_rs_rate2[i]->Sumw2();
110  h_clawss_rs_rate1W[i]->Sumw2();
111  h_clawss_rs_rate2W[i]->Sumw2();
112  h_clawss_pe1[i]->Sumw2();
113  h_clawss_pe2[i]->Sumw2();
114  h_clawss_pe1W[i]->Sumw2();
115  h_clawss_pe2W[i]->Sumw2();
116  }
117 }
118 
119 
120 void ClawsStudyModule::initialize()
121 {
122  B2INFO("ClawsStudyModule: Initialize");
123 
124  REG_HISTOGRAM
125 
126  //get CLAWSS paramters ie energy threshold
127  getXMLData();
128 }
129 
130 void ClawsStudyModule::beginRun()
131 {
132 }
133 
134 void ClawsStudyModule::event()
135 {
136  //Here comes the actual event processing
137  StoreArray<CLAWSSimHit> SimHits;
139  StoreArray<SADMetaHit> MetaHits;
140 
141  double rate = 0;
142  int ring_section = -1;
143  for (const auto& MetaHit : MetaHits) {
144  rate = MetaHit.getrate();
145  ring_section = MetaHit.getring_section() - 1;
146  }
147 
148  /*
149  //number of entries in SimHits
150  int nSimHits = SimHits.getEntries();
151 
152  //loop over all SimHit entries
153  for (int i = 0; i < nSimHits; i++) {
154  CLAWSSimHit* aHit = SimHits[i];
155  int lad = aHit->getLadder();
156  int sen = aHit->getSensor();
157  //const int detNb = SimHit.getCellId();
158  //int pdg = SimHit.getPDGCode();
159  int detNB = (lad - 1) * 8 + sen - 1;
160  //int detNB = aHit->getCellId();
161  if (detNB < 16) {
162  //int trkID = aHit->getTrackId();
163  int pdg = aHit->getPDG();
164  double Edep = aHit->getEnergyVisible() * 1e3; //GeV -> MeV
165  double tof = aHit->getTime(); //ns
166 
167  //h_clawss_Evtof1[detNB]->Fill(tof, Edep);
168  //if (pdg == 22) h_clawss_Evtof2[detNB]->Fill(tof, Edep);
169  //else if (fabs(pdg) == 11) h_clawss_Evtof3[detNB]->Fill(tof, Edep);
170  //else h_clawss_Evtof4[detNB]->Fill(tof, Edep);
171  if (Edep > m_Ethres) {
172  //h_clawss_edep[detNB]->Fill(Edep);
173  //h_Wclawss_edep[detNB]->Fill(Edep, rate);
174  }
175 
176  }
177  }
178  */
179  for (const auto& Hit : Hits) {
180  const int detNb = Hit.getdetNb();
181  if (detNb < 16) {
182  const int timebin = Hit.gettime();
183  const float edep = Hit.getedep();
184  const float pe = Hit.getPE();
185  h_clawss_hitrate1->Fill(detNb);
186  h_clawss_hitrate1W->Fill(detNb, rate);
187  h_clawss_rate1[detNb]->Fill(pe);
188  h_clawss_rate1W[detNb]->Fill(pe, rate);
189  h_clawss_rs_rate1[detNb]->Fill(pe, ring_section);
190  h_clawss_rs_rate1W[detNb]->Fill(pe, ring_section, rate);
191  h_clawss_rs_hitrate1->Fill(detNb, ring_section);
192  h_clawss_rs_hitrate1W->Fill(detNb, ring_section, rate);
193  h_clawss_pe1[detNb]->Fill(timebin, pe);
194  h_clawss_pe1W[detNb]->Fill(timebin, pe, rate);
195  if (edep > m_Ethres) {
196  h_clawss_hitrate2->Fill(detNb);
197  h_clawss_hitrate2W->Fill(detNb, rate);
198  h_clawss_rate2[detNb]->Fill(pe);
199  h_clawss_rate2W[detNb]->Fill(pe, rate);
200  h_clawss_rs_rate2[detNb]->Fill(pe, ring_section);
201  h_clawss_rs_rate2W[detNb]->Fill(pe, ring_section, rate);
202  h_clawss_rs_hitrate2->Fill(detNb, ring_section);
203  h_clawss_rs_hitrate2W->Fill(detNb, ring_section, rate);
204  h_clawss_pe2[detNb]->Fill(timebin, pe);
205  h_clawss_pe2W[detNb]->Fill(timebin, pe, rate);
206  }
207  }
208  }
209 
210 }
211 
212 //read energy threshold from CLAWS.xml
213 void ClawsStudyModule::getXMLData()
214 {
215  GearDir content = GearDir("/Detector/DetectorComponent[@name=\"CLAWS\"]/Content/");
216  m_Ethres = content.getDouble("Ethres");
217 
218  B2INFO("ClawsStudy");
219 }
220 
221 
222 void ClawsStudyModule::endRun()
223 {
224 
225 
226 
227 }
228 
229 void ClawsStudyModule::terminate()
230 {
231 }
232 
233 
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::claws::ClawsStudyModule
Study module for Clawss (BEAST)
Definition: ClawsStudyModule.h:39
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::GearDir
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:41
Belle2::StoreArray
Accessor to arrays stored in the data store.
Definition: ECLMatchingPerformanceExpertModule.h:33
Belle2::HistoModule
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
Belle2::Hit
Structure to hold some of the calpulse data.
Definition: TOPTimeBaseCalibratorModule.h:40