Belle II Software development
ClawStudyModule.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#include <beast/claw/modules/ClawStudyModule.h>
10#include <beast/claw/dataobjects/ClawSimHit.h>
11#include <beast/claw/dataobjects/ClawHit.h>
12#include <generators/SAD/dataobjects/SADMetaHit.h>
13#include <framework/datastore/StoreArray.h>
14#include <framework/gearbox/GearDir.h>
15#include <framework/gearbox/Const.h>
16#include <framework/logging/Logger.h>
17#include <cmath>
18
19#include <fstream>
20#include <string>
21
22// ROOT
23#include <TH1.h>
24#include <TH2.h>
25
26int eventNum = 0;
27
28using namespace std;
29
30using namespace Belle2;
31using namespace claw;
32
33//-----------------------------------------------------------------
34// Register the Module
35//-----------------------------------------------------------------
36REG_MODULE(ClawStudy);
37
38//-----------------------------------------------------------------
39// Implementation
40//-----------------------------------------------------------------
41
43{
44 // Set module properties
45 setDescription("Study module for Claws (BEAST)");
46
47 addParam("Ethres", m_Ethres, "Energy threshold in MeV", 0.0);
48}
49
51{
52}
53
54//This module is a histomodule. Any histogram created here will be saved by the HistoManager module
56{
57 for (int i = 0; i < 8; i++) {
58 h_claws_Evtof1[i] = new TH2F(TString::Format("claws_Evtof1_%d", i), "Energy deposited [MeV] vs TOF [ns] - all", 5000, 0., 1000.,
59 1000, 0., 10.);
60 h_claws_Evtof2[i] = new TH2F(TString::Format("claws_Evtof2_%d", i), "Energy deposited [MeV] vs TOF [ns] - only photons", 5000, 0.,
61 1000., 1000, 0., 10.);
62 h_claws_Evtof3[i] = new TH2F(TString::Format("claws_Evtof3_%d", i), "Energy deposited [MeV] vs TOF [ns] - only e+/e-", 5000, 0.,
63 1000., 1000, 0., 10.);
64 h_claws_Evtof4[i] = new TH2F(TString::Format("claws_Evtof4_%d", i), "Energy deposited [MeV] vs TOF [ns] - only e+/e-", 5000, 0.,
65 1000., 1000, 0., 10.);
66 h_claws_edep[i] = new TH1F(TString::Format("claws_edep_%d", i), "Energy deposited [MeV]", 5000, 0., 10.);
67 h_Wclaws_edep[i] = new TH1F(TString::Format("Wclaws_edep_%d", i), "Energy deposited [MeV]", 5000, 0., 10.);
68 }
69
70 h_claws_hitrate1 = new TH1F("claws_hitrate1", "Hit distributions", 8, 0., 8.);
71 h_claws_hitrate2 = new TH1F("claws_hitrate2", "Hit distributions", 8, 0., 8.);
72 h_claws_hitrate1W = new TH1F("claws_hitrate1W", "Hit distributions", 8, 0., 8.);
73 h_claws_hitrate2W = new TH1F("claws_hitrate2W", "Hit distributions", 8, 0., 8.);
74
75 h_claws_hitrate1->Sumw2();
76 h_claws_hitrate1W->Sumw2();
77 h_claws_hitrate2->Sumw2();
78 h_claws_hitrate2W->Sumw2();
79
80 h_claws_rs_hitrate1 = new TH2F("claws_rs_hitrate1", "Hit distributions vs rs", 8, 0., 8., 12, 0., 12.);
81 h_claws_rs_hitrate2 = new TH2F("claws_rs_hitrate2", "Hit distributions vs rs", 8, 0., 8., 12, 0., 12.);
82 h_claws_rs_hitrate1W = new TH2F("claws_rs_hitrate1W", "Hit distributions vs rs", 8, 0., 8., 12, 0., 12.);
83 h_claws_rs_hitrate2W = new TH2F("claws_rs_hitrate2W", "Hit distributions vs rs", 8, 0., 8., 12, 0., 12.);
84
85 h_claws_rs_hitrate1->Sumw2();
86 h_claws_rs_hitrate1W->Sumw2();
87 h_claws_rs_hitrate2->Sumw2();
88 h_claws_rs_hitrate2W->Sumw2();
89
90 for (int i = 0; i < 8; i++) {
91 h_claws_rate1[i] = new TH1F(TString::Format("claws_rate1_%d", i), "PE distributions", 5000, 0., 5000.);
92 h_claws_rate2[i] = new TH1F(TString::Format("claws_rate2_%d", i), "PE distributions", 5000, 0., 5000.);
93 h_claws_rate1W[i] = new TH1F(TString::Format("claws_rate1W_%d", i), "PE distributions", 5000, 0., 5000.);
94 h_claws_rate2W[i] = new TH1F(TString::Format("claws_rate2W_%d", i), "PE distributions", 5000, 0., 5000.);
95 h_claws_pe1[i] = new TH2F(TString::Format("claws_pe1_%d", i), "PE distributions", 5000, 0., 5000., 1000, 0., 1000.);
96 h_claws_pe2[i] = new TH2F(TString::Format("claws_pe2_%d", i), "PE distributions", 5000, 0., 5000., 1000, 0., 1000.);
97 h_claws_pe1W[i] = new TH2F(TString::Format("claws_pe1W_%d", i), "PE distributions", 5000, 0., 5000., 1000, 0., 1000.);
98 h_claws_pe2W[i] = new TH2F(TString::Format("claws_pe2W_%d", i), "PE distributions", 5000, 0., 5000., 1000, 0., 1000.);
99
100 h_claws_rs_rate1[i] = new TH2F(TString::Format("claws_rs_rate1_%d", i), "PE distributions", 5000, 0., 5000., 12, 0., 12.);
101 h_claws_rs_rate2[i] = new TH2F(TString::Format("claws_rs_rate2_%d", i), "PE distributions", 5000, 0., 5000., 12, 0., 12.);
102 h_claws_rs_rate1W[i] = new TH2F(TString::Format("claws_rs_rate1W_%d", i), "PE distributions", 5000, 0., 5000., 12, 0., 12.);
103 h_claws_rs_rate2W[i] = new TH2F(TString::Format("claws_rs_rate2W_%d", i), "PE distributions", 5000, 0., 5000., 12, 0., 12.);
104
105 h_claws_rate1[i]->Sumw2();
106 h_claws_rate2[i]->Sumw2();
107 h_claws_rate1W[i]->Sumw2();
108 h_claws_rate2W[i]->Sumw2();
109 h_claws_rs_rate1[i]->Sumw2();
110 h_claws_rs_rate2[i]->Sumw2();
111 h_claws_rs_rate1W[i]->Sumw2();
112 h_claws_rs_rate2W[i]->Sumw2();
113 h_claws_pe1[i]->Sumw2();
114 h_claws_pe2[i]->Sumw2();
115 h_claws_pe1W[i]->Sumw2();
116 h_claws_pe2W[i]->Sumw2();
117 }
118}
119
120
122{
123 B2INFO("ClawStudyModule: Initialize");
124
125 REG_HISTOGRAM
126
127 //get CLAWS paramters ie energy threshold
128 getXMLData();
129}
130
132{
133}
134
136{
137 //Here comes the actual event processing
140 StoreArray<SADMetaHit> MetaHits;
141
142 double rate = 0;
143 int ring_section = -1;
144 for (const auto& MetaHit : MetaHits) {
145 rate = MetaHit.getrate();
146 ring_section = MetaHit.getring_section() - 1;
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 < (int) SimHits.getEntries(); i++) {
154 ClawSimHit* aHit = SimHits[i];
155 int detNB = aHit->getCellId();
156 if (detNB < 8) {
157 //int trkID = aHit->getTrackId();
158 int pdg = aHit->getPDGCode();
159 double Edep = aHit->getEnergyDep() * 1e3; //GeV -> MeV
160 double tof = aHit->getFlightTime(); //ns
161
162 h_claws_Evtof1[detNB]->Fill(tof, Edep);
163 if (pdg == Const::photon.getPDGCode()) h_claws_Evtof2[detNB]->Fill(tof, Edep);
164 else if (fabs(pdg) == Const::electron.getPDGCode()) h_claws_Evtof3[detNB]->Fill(tof, Edep);
165 else h_claws_Evtof4[detNB]->Fill(tof, Edep);
166 if (Edep > m_Ethres) {
167 h_claws_edep[detNB]->Fill(Edep);
168 h_Wclaws_edep[detNB]->Fill(Edep, rate);
169 }
170 }
171 }
172
173 for (const auto& Hit : Hits) {
174 const int detNb = Hit.getdetNb();
175 if (detNb < 8) {
176 const int timebin = Hit.gettime();
177 const float edep = Hit.getedep();
178 const float pe = Hit.getPE();
179 h_claws_hitrate1->Fill(detNb);
180 h_claws_hitrate1W->Fill(detNb, rate);
181 h_claws_rate1[detNb]->Fill(pe);
182 h_claws_rate1W[detNb]->Fill(pe, rate);
183 h_claws_rs_rate1[detNb]->Fill(pe, ring_section);
184 h_claws_rs_rate1W[detNb]->Fill(pe, ring_section, rate);
185 h_claws_rs_hitrate1->Fill(detNb, ring_section);
186 h_claws_rs_hitrate1W->Fill(detNb, ring_section, rate);
187 h_claws_pe1[detNb]->Fill(timebin, pe);
188 h_claws_pe1W[detNb]->Fill(timebin, pe, rate);
189 if (edep > m_Ethres) {
190 h_claws_hitrate2->Fill(detNb);
191 h_claws_hitrate2W->Fill(detNb, rate);
192 h_claws_rate2[detNb]->Fill(pe);
193 h_claws_rate2W[detNb]->Fill(pe, rate);
194 h_claws_rs_rate2[detNb]->Fill(pe, ring_section);
195 h_claws_rs_rate2W[detNb]->Fill(pe, ring_section, rate);
196 h_claws_rs_hitrate2->Fill(detNb, ring_section);
197 h_claws_rs_hitrate2W->Fill(detNb, ring_section, rate);
198 h_claws_pe2[detNb]->Fill(timebin, pe);
199 h_claws_pe2W[detNb]->Fill(timebin, pe, rate);
200 }
201 }
202 }
203
204}
205
206//read energy threshold from CLAW.xml
208{
209 GearDir content = GearDir("/Detector/DetectorComponent[@name=\"CLAW\"]/Content/");
210 m_Ethres = content.getDouble("Ethres");
211
212 B2INFO("ClawStudy");
213}
214
215
217{
218
219
220
221}
222
224{
225}
226
227
ClassClawSimHit - Geant4 simulated hit for the Claw crystal in beast.
Definition: ClawSimHit.h:31
int getPDGCode() const
Get Particle PDG (can be one of secondaries)
Definition: ClawSimHit.h:99
int getCellId() const
Get Cell ID.
Definition: ClawSimHit.h:89
double getFlightTime() const
Get Flight time from IP.
Definition: ClawSimHit.h:104
double getEnergyDep() const
Get Deposit energy.
Definition: ClawSimHit.h:109
int getPDGCode() const
PDG code.
Definition: Const.h:473
static const ParticleType photon
photon particle
Definition: Const.h:673
static const ChargedStable electron
electron particle
Definition: Const.h:659
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:31
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
TH2F * h_claws_pe2[8]
Energy deposited.
virtual ~ClawStudyModule()
Destructor.
TH2F * h_claws_rs_hitrate1
Energy deposited.
TH2F * h_claws_rs_hitrate2
Energy deposited.
ClawStudyModule()
Constructor: Sets the description, the properties and the parameters of the module.
TH2F * h_claws_rs_rate2W[8]
Energy deposited.
TH1F * h_claws_rate1[8]
Energy deposited.
virtual void initialize() override
Initialize the Module.
TH2F * h_claws_pe1W[8]
Energy deposited.
TH1F * h_claws_hitrate1
Energy deposited.
TH2F * h_claws_rs_hitrate1W
Energy deposited.
virtual void event() override
Event processor.
TH2F * h_claws_rs_hitrate2W
Energy deposited.
TH1F * h_claws_hitrate2
Energy deposited.
TH2F * h_claws_rs_rate1W[8]
Energy deposited.
TH1F * h_claws_edep[8]
Energy deposited.
virtual void endRun() override
End-of-run action.
TH1F * h_claws_hitrate2W
Energy deposited.
virtual void getXMLData()
reads data from CLAW.xml
TH2F * h_claws_rs_rate2[8]
Energy deposited.
TH2F * h_claws_Evtof4[8]
Energy deposited vs TOF.
virtual void terminate() override
Termination action.
virtual void beginRun() override
Called when entering a new run.
TH1F * h_claws_rate2[8]
Energy deposited.
TH2F * h_claws_Evtof3[8]
Energy deposited vs TOF.
TH2F * h_claws_pe2W[8]
Energy deposited.
TH2F * h_claws_rs_rate1[8]
Energy deposited.
TH1F * h_Wclaws_edep[8]
Energy deposited.
TH2F * h_claws_Evtof2[8]
Energy deposited vs TOF.
TH1F * h_claws_hitrate1W
Energy deposited.
TH1F * h_claws_rate2W[8]
Energy deposited.
TH2F * h_claws_Evtof1[8]
Energy deposited vs TOF.
double m_Ethres
Energy threshold.
TH1F * h_claws_rate1W[8]
Energy deposited.
TH2F * h_claws_pe1[8]
Energy deposited.
virtual void defineHisto() override
Defines the histograms.
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.
STL namespace.
Structure to hold some of the calpulse data.