Belle II Software development
FANGSStudyModule.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/fangs/modules/FANGSStudyModule.h>
10#include <beast/fangs/dataobjects/FANGSSimHit.h>
11#include <beast/fangs/dataobjects/FANGSHit.h>
12#include <framework/datastore/StoreArray.h>
13#include <framework/gearbox/GearDir.h>
14#include <framework/logging/Logger.h>
15#include <framework/gearbox/Const.h>
16#include <cmath>
17
18#include <string>
19
20// ROOT
21#include <Math/Vector3D.h>
22
23int eventNum = 0;
24
25using namespace std;
26
27using namespace Belle2;
28using namespace fangs;
29
30//-----------------------------------------------------------------
31// Register the Module
32//-----------------------------------------------------------------
33REG_MODULE(FANGSStudy);
34
35//-----------------------------------------------------------------
36// Implementation
37//-----------------------------------------------------------------
38
40{
41 // Set module properties
42 setDescription("Study module for Fangs (BEAST)");
43
44}
45
47{
48}
49
50//This module is a histomodule. Any histogram created here will be saved by the HistoManager module
52{
53 h_time = new TH2F("h_time", "Detector # vs. time", 20, 0., 20., 1000, 0., 750.);
54 h_time->Sumw2();
55 h_timeWeighted = new TH2F("h_timeWeigthed", "Detector # vs. time weighted by the energy deposited", 20, 0., 20., 1000, 0., 750.);
56 h_timeWeighted->Sumw2();
57 h_timeThres = new TH2F("h_timeThres", "Detector # vs. time", 20, 0., 20., 750, 0., 750.);
58 h_timeThres->Sumw2();
59 h_timeWeightedThres = new TH2F("h_timeWeigthedThres", "Detector # vs. time weighted by the energy deposited", 20, 0., 20., 750, 0.,
60 750.);
61 h_timeWeightedThres->Sumw2();
62 h_edep = new TH2F("h_edep", "Time bin # vs. energy deposited", 20, 0., 20., 1000, 0., 10.);
63 h_edep->Sumw2();
64 h_edep1 = new TH2F("h_edep1", "Time bin # vs. energy deposited", 20, 0., 20., 1000, 0., 10.);
65 h_edep1->Sumw2();
66 h_edep2 = new TH2F("h_edep2", "Time bin # vs. energy deposited", 20, 0., 20., 1000, 0., 10.);
67 h_edep2->Sumw2();
68 h_edep3 = new TH2F("h_edep3", "Time bin # vs. energy deposited", 20, 0., 20., 1000, 0., 10.);
69 h_edep3->Sumw2();
70
71 h_edepThres = new TH2F("h_edepThres", "Time bin # vs. energy deposited", 20, 0., 20., 1000, 0., 10.);
72 h_edepThres->Sumw2();
73 h_edepThres1 = new TH2F("h_edepThres1", "Time bin # vs. energy deposited", 20, 0., 20., 1000, 0., 10.);
74 h_edepThres1->Sumw2();
75 h_edepThres2 = new TH2F("h_edepThres2", "Time bin # vs. energy deposited", 20, 0., 20., 1000, 0., 10.);
76 h_edepThres2->Sumw2();
77 h_edepThres3 = new TH2F("h_edepThres3", "Time bin # vs. energy deposited", 20, 0., 20., 1000, 0., 10.);
78 h_edepThres3->Sumw2();
79 for (int i = 0; i < 3; i++) {
80 h_zvedep[i] = new TH1F(TString::Format("h_zvedep_%d", i), "edep [MeV] vs. z [cm]", 2000, -25., 25.);
81 h_zvedep[i]->Sumw2();
82
83 h_xvzvedep[i] = new TH2F(TString::Format("h_xvzvedep_%d", i), "edep [MeV] vs. z [cm]", 2000, -25., 25., 2000, -25., 25.);
84 h_xvzvedep[i]->Sumw2();
85
86 h_yvzvedep[i] = new TH2F(TString::Format("h_yvzvedep_%d", i), "edep [MeV] vs. z [cm]", 2000, -25., 25., 2000, -25., 25.);
87 h_yvzvedep[i]->Sumw2();
88
89 h_rvzvedep[i] = new TH2F(TString::Format("h_rvzvedep_%d", i), "edep [MeV] vs. z [cm]", 2000, 0., 25., 2000, -25., 25.);
90 h_rvzvedep[i]->Sumw2();
91
92 h_xvyvedep[i] = new TH2F(TString::Format("h_xvyvedep_%d", i), "edep [MeV] vs. y [cm]", 2000, -25., 25., 2000, -25., 25.);
93 h_xvyvedep[i]->Sumw2();
94
95 h_zvedepW[i] = new TH1F(TString::Format("h_zvedepW_%d", i), "edep [MeV] vs. z [cm]", 2000, -25., 25.);
96 h_zvedepW[i]->Sumw2();
97
98 h_xvzvedepW[i] = new TH2F(TString::Format("h_xvzvedepW_%d", i), "edep [MeV] vs. x vs. z [cm]", 2000, -25., 25., 2000, -25., 25.);
99 h_xvzvedepW[i]->Sumw2();
100
101 h_yvzvedepW[i] = new TH2F(TString::Format("h_yvzvedepW_%d", i), "edep [MeV] vs. y vs. z [cm]", 2000, -25., 25., 2000, -25., 25.);
102 h_yvzvedepW[i]->Sumw2();
103
104 h_xvyvedepW[i] = new TH2F(TString::Format("h_xvyvedepW_%d", i), "edep [MeV] vs. x vs. y [cm]", 2000, -25., 25., 2000, -25., 25.);
105 h_xvyvedepW[i]->Sumw2();
106
107 h_rvzvedepW[i] = new TH2F(TString::Format("h_rvzvedepW_%d", i), "edep [MeV] vs. z [cm]", 2000, 0., 25., 2000, -25., 25.);
108 h_rvzvedepW[i]->Sumw2();
109
110
111 h_zvedepT[i] = new TH1F(TString::Format("h_zvedepT_%d", i), "edep [MeV] vs. z [cm]", 2000, -25., 25.);
112 h_zvedepT[i]->Sumw2();
113
114 h_xvzvedepT[i] = new TH2F(TString::Format("h_xvzvedepT_%d", i), "edep [MeV] vs. z [cm]", 2000, -25., 25., 2000, -25., 25.);
115 h_xvzvedepT[i]->Sumw2();
116
117 h_yvzvedepT[i] = new TH2F(TString::Format("h_yvzvedepT_%d", i), "edep [MeV] vs. z [cm]", 2000, -25., 25., 2000, -25., 25.);
118 h_yvzvedepT[i]->Sumw2();
119
120 h_rvzvedepT[i] = new TH2F(TString::Format("h_rvzvedepT_%d", i), "edep [MeV] vs. z [cm]", 2000, 0., 25., 2000, -25., 25.);
121 h_rvzvedepT[i]->Sumw2();
122
123 h_xvyvedepT[i] = new TH2F(TString::Format("h_xvyvedepT_%d", i), "edep [MeV] vs. y [cm]", 2000, -25., 25., 2000, -25., 25.);
124 h_xvyvedepT[i]->Sumw2();
125
126 h_zvedepWT[i] = new TH1F(TString::Format("h_zvedepWT_%d", i), "edep [MeV] vs. z [cm]", 2000, -25., 25.);
127 h_zvedepWT[i]->Sumw2();
128
129 h_xvzvedepWT[i] = new TH2F(TString::Format("h_xvzvedepWT_%d", i), "edep [MeV] vs. x vs. z [cm]", 2000, -25., 25., 2000, -25., 25.);
130 h_xvzvedepWT[i]->Sumw2();
131
132 h_yvzvedepWT[i] = new TH2F(TString::Format("h_yvzvedepWT_%d", i), "edep [MeV] vs. y vs. z [cm]", 2000, -25., 25., 2000, -25., 25.);
133 h_yvzvedepWT[i]->Sumw2();
134
135 h_xvyvedepWT[i] = new TH2F(TString::Format("h_xvyvedepWT_%d", i), "edep [MeV] vs. x vs. y [cm]", 2000, -25., 25., 2000, -25., 25.);
136 h_xvyvedepWT[i]->Sumw2();
137
138 h_rvzvedepWT[i] = new TH2F(TString::Format("h_rvzvedepWT_%d", i), "edep [MeV] vs. z [cm]", 2000, 0., 25., 2000, -25., 25.);
139 h_rvzvedepWT[i]->Sumw2();
140 }
141 h_Edep = new TH2F("h_Edep", "det # # vs. energy deposited", 20, 0., 20., 1000, 0., 10.);
142 h_pxNb = new TH2F("h_pxNb", "det # # vs. nb pixel", 20, 0., 20., 1000, 0., 1000.);
143 for (int i = 0; i < 15; i++) {
144 h_cvr[i] = new TH2F(TString::Format("cvr_%d", i), " col v. row", 80, 0., 80., 336, 0., 336.);
145 }
146}
147
148
150{
151 B2INFO("FANGSStudyModule: Initialize");
152
153 REG_HISTOGRAM
154
155 //convert sample time into rate in Hz
156 //rateCorrection = m_sampletime / 1e6;
157 //get FANGS parameters
158 getXMLData();
159
160 fctQ_Calib1 = new TF1("fctQ_Calib1", "[0]*([1]*x-[2])/([3]-x)", 0., 15.);
161 fctQ_Calib1->SetParameters(m_TOTQ1, m_TOTC1, m_TOTA1 * m_TOTB1, m_TOTA1);
162
163 fctQ_Calib2 = new TF1("fctQ_Calib2", "[0]*([1]*x-[2])/([3]-x)", 0., 15.);
164 fctQ_Calib2->SetParameters(m_TOTQ2, m_TOTC2, m_TOTA2 * m_TOTB2, m_TOTA2);
165}
166
168{
169}
170
172{
173 //Here comes the actual event processing
174
177
178 int olddetNb = -1;
179 int ipix = 0;
180 float esum = 0;
181 //number of entries in Hits
182 for (const auto& FANGSHit : Hits) {
183 int detNb = FANGSHit.getdetNb();
184 //int pdg = FANGSHit.getPDG();
185 //int trkID = FANGSHit.gettrkID();
186 int col = FANGSHit.getcolumn();
187 int row = FANGSHit.getrow();
188 int tot = FANGSHit.getTOT();
189 int bcid = FANGSHit.getBCID();
190
191 if (olddetNb != detNb) {
192 if (esum > 0) {
193 h_Edep->Fill(detNb, esum);
194 h_pxNb->Fill(detNb, ipix);
195 }
196 ipix = 0;
197 esum = 0;
198 for (int j = 0; j < maxSIZE; j++) {
199 x[j] = 0;
200 y[j] = 0;
201 z[j] = 0;
202 e[j] = 0;
203 }
204 olddetNb = detNb;
205 }
206 x[ipix] = col * (2. * m_ChipColumnX / (float)m_ChipColumnNb) - m_ChipColumnX;
207 y[ipix] = row * (2. * m_ChipRowY / (float)m_ChipRowNb) - m_ChipRowY;
208 z[ipix] = (m_PixelTimeBin / 2. + m_PixelTimeBin * bcid) * m_v_sensor;
209 if (tot < 3) e[ipix] = fctQ_Calib1->Eval(tot) * m_Workfct * 1e-3;
210 else e[ipix] = fctQ_Calib2->Eval(tot) * m_Workfct * 1e-3;
211 esum += e[ipix];
212 h_cvr[detNb]->Fill(col, row);
213 ipix ++;
214 }
215 //number of entries in SimHits
216 int nSimHits = SimHits.getEntries();
217 //cout << nSimHits << endl;
218
219 //loop over all SimHit entries
220 for (int i = 0; i < nSimHits; i++) {
221 FANGSSimHit* aHit = SimHits[i];
222 int lad = aHit->getLadder();
223 int sen = aHit->getSensor();
224 double adep = aHit->getEnergyDep();
225 double timeBin = aHit->getTime();
226 int pdg = aHit->getPDG();
227
228 ROOT::Math::XYZVector position = aHit->getPosEntry();
229 double r = sqrt(position.X() * position.X() + position.Y() * position.Y());
230 int detNB = (lad - 1) * 5 + sen - 1;
231 //cout <<" lad " << lad << " sen " << sen << " detNB " << detNB << " time " << timeBin << " edep " << adep*1e3 << endl;
232 //cout <<" x " << position.X() << " y " << position.Y() << " z " << position.Z() << endl;
233 h_time->Fill(detNB, timeBin);
234 h_edep->Fill(detNB, adep * 1e3);
235 if (fabs(pdg) == Const::electron.getPDGCode())h_edep1->Fill(detNB, adep * 1e3);
236 if (pdg == Const::photon.getPDGCode())h_edep2->Fill(detNB, adep * 1e3);
237 if (pdg != Const::photon.getPDGCode() && fabs(pdg) != Const::electron.getPDGCode())h_edep3->Fill(detNB, adep * 1e3);
238 if (adep > 50.*1e-6) {
239 h_timeThres->Fill(detNB, timeBin);
240 h_edepThres->Fill(detNB, adep * 1e3);
241 if (fabs(pdg) == Const::electron.getPDGCode())h_edepThres1->Fill(detNB, adep * 1e3);
242 if (pdg == Const::photon.getPDGCode())h_edepThres2->Fill(detNB, adep * 1e3);
243 if (pdg != Const::photon.getPDGCode() && fabs(pdg) != Const::electron.getPDGCode())h_edepThres3->Fill(detNB, adep * 1e3);
244 }
245 h_zvedep[lad - 1]->Fill(position.Z());
246 h_xvzvedep[lad - 1]->Fill(position.X(), position.Z());
247 h_yvzvedep[lad - 1]->Fill(position.Y(), position.Z());
248 h_xvyvedep[lad - 1]->Fill(position.X(), position.Y());
249 h_rvzvedep[lad - 1]->Fill(r, position.Z());
250 h_zvedepW[lad - 1]->Fill(position.Z(), adep * 1e3);
251 h_xvzvedepW[lad - 1]->Fill(position.X(), position.Z(), adep * 1e3);
252 h_yvzvedepW[lad - 1]->Fill(position.Y(), position.Z(), adep * 1e3);
253 h_xvyvedepW[lad - 1]->Fill(position.X(), position.Y(), adep * 1e3);
254 h_rvzvedepW[lad - 1]->Fill(r, position.Z(), adep * 1e3);
255 if (adep > 50.*1e-6) {
256 h_zvedepT[lad - 1]->Fill(position.Z());
257 h_xvzvedepT[lad - 1]->Fill(position.X(), position.Z());
258 h_yvzvedepT[lad - 1]->Fill(position.Y(), position.Z());
259 h_xvyvedepT[lad - 1]->Fill(position.X(), position.Y());
260 h_rvzvedepT[lad - 1]->Fill(r, position.Z());
261 h_zvedepWT[lad - 1]->Fill(position.Z(), adep * 1e3);
262 h_xvzvedepWT[lad - 1]->Fill(position.X(), position.Z(), adep * 1e3);
263 h_yvzvedepWT[lad - 1]->Fill(position.Y(), position.Z(), adep * 1e3);
264 h_xvyvedepWT[lad - 1]->Fill(position.X(), position.Y(), adep * 1e3);
265 h_rvzvedepWT[lad - 1]->Fill(r, position.Z(), adep * 1e3);
266 }
267 }
268
269 eventNum++;
270}
271
273{
274 GearDir content = GearDir("/Detector/DetectorComponent[@name=\"FANGS\"]/Content/");
275
276 m_PixelThreshold = content.getInt("PixelThreshold");
277 m_PixelThresholdRMS = content.getInt("PixelThresholdRMS");
278 m_PixelTimeBinNb = content.getInt("PixelTimeBinNb");
279 m_PixelTimeBin = content.getDouble("PixelTimeBin");
280 m_ChipColumnNb = content.getInt("ChipColumnNb");
281 m_ChipRowNb = content.getInt("ChipRowNb");
282 m_ChipColumnX = content.getDouble("ChipColumnX");
283 m_ChipRowY = content.getDouble("ChipRowY");
284 m_TOTA1 = content.getDouble("TOTA1");
285 m_TOTB1 = content.getDouble("TOTB1");
286 m_TOTC1 = content.getDouble("TOTC1");
287 m_TOTQ1 = content.getDouble("TOTQ1");
288 m_TOTA2 = content.getDouble("TOTA2");
289 m_TOTB2 = content.getDouble("TOTB2");
290 m_TOTC2 = content.getDouble("TOTC2");
291 m_TOTQ2 = content.getDouble("TOTQ2");
292 m_v_sensor = content.getDouble("v_sensor");
293 m_sensor_width = content.getDouble("sensor_width");
294 m_Workfct = content.getDouble("Workfct");
295
296}
297
299{
300
301
302
303}
304
306{
307}
308
309
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
ClassFANGSHit - digitization simulated hit for the FANGS detector.
Definition: FANGSHit.h:26
int getBCID() const
Return the BCID.
Definition: FANGSHit.h:54
int getrow() const
Return the row.
Definition: FANGSHit.h:52
int getTOT() const
Return the TOT.
Definition: FANGSHit.h:56
int getdetNb() const
Return the TPC number.
Definition: FANGSHit.h:58
int getcolumn() const
Return the column.
Definition: FANGSHit.h:50
Class FANGSSimHit - Geant4 simulated hit for the FANGS detector.
Definition: FANGSSimHit.h:28
ROOT::Math::XYZVector getPosEntry() const
Return the entry track position.
Definition: FANGSSimHit.h:78
float getTime() const
Return the global time.
Definition: FANGSSimHit.h:74
float getEnergyDep() const
Return the energy deposition in electrons.
Definition: FANGSSimHit.h:76
int getLadder() const
Return the Ladder number (starting at 1, increasing with phi)
Definition: FANGSSimHit.h:68
int getSensor() const
Return the Sensor number (starting at 1, increasing with decreasing z)
Definition: FANGSSimHit.h:70
int getPDG() const
Return the PDG number of the track.
Definition: FANGSSimHit.h:72
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_timeThres
Time distribution with energy threshold applied.
double m_ChipRowY
Chip row y dimension.
TH2F * h_rvzvedepW[3]
Energy vs x vs y.
TH2F * h_yvzvedepWT[3]
Energy vs y vs z.
TH2F * h_edepThres1
Energy deposited above threshold per time bin.
TH2F * h_edep
Energy deposited per time bin.
TH2F * h_edep1
Energy deposited per time bin.
double m_PixelTimeBin
Pixel time bin.
TH2F * h_timeWeighted
Time distribution weighted per the energy deposited.
virtual void initialize() override
Initialize the Module.
virtual void event() override
Event processor.
TH2F * h_edep3
Energy deposited per time bin.
TH2F * h_time
Time distribution.
virtual void endRun() override
End-of-run action.
double m_ChipColumnX
Chip column x dimension.
virtual void getXMLData()
reads data from FANGS.xml
virtual void terminate() override
Termination action.
double m_v_sensor
Drift velocity in sensor.
TH2F * h_yvzvedepT[3]
Energy vs y vs z.
TH2F * h_xvzvedepT[3]
Energy vs x vs z.
FANGSStudyModule()
Constructor: Sets the description, the properties and the parameters of the module.
TH2F * h_yvzvedep[3]
Energy vs y vs z.
TH2F * h_edep2
Energy deposited per time bin.
TH1F * h_zvedepW[3]
Energy deposited vs z.
TH2F * h_xvzvedep[3]
Energy vs x vs z.
TH2F * h_rvzvedepT[3]
Energy vs x vs y.
TH2F * h_xvyvedepWT[3]
Energy vs x vs y.
TH2F * h_Edep
Digitized energy deposited per detector.
TH1F * h_zvedepWT[3]
Energy deposited vs z.
float y[maxSIZE]
y point of the track
TH1F * h_zvedep[3]
Energy deposited vs z.
TH1F * h_zvedepT[3]
Energy deposited vs z.
virtual ~FANGSStudyModule()
Destructor.
virtual void beginRun() override
Called when entering a new run.
TH2F * h_rvzvedepWT[3]
Energy vs x vs y.
TH2F * h_timeWeightedThres
Time distribution weighted per the energy deposited with energy threshold applied.
TH2F * h_edepThres2
Energy deposited above threshold per time bin.
TH2F * h_edepThres
Energy deposited above threshold per time bin.
TH2F * h_xvyvedep[3]
Energy vs x vs y.
TH2F * h_rvzvedep[3]
Energy vs x vs y.
TH2F * h_yvzvedepW[3]
Energy vs y vs z.
TF1 * fctQ_Calib2
Define Q calib 2.
TH2F * h_xvzvedepW[3]
Energy vs x vs z.
int m_PixelTimeBinNb
Pixel time number of bin.
int m_ChipColumnNb
Chip column number.
TH2F * h_xvzvedepWT[3]
Energy vs x vs z.
TH2F * h_pxNb
Pixel number per detector.
TF1 * fctQ_Calib1
Define Q calib 1.
TH2F * h_xvyvedepW[3]
Energy vs x vs y.
int m_PixelThresholdRMS
Pixel threshold RMS.
int m_PixelThreshold
Pixel threshold.
float e[maxSIZE]
e point of the track
float x[maxSIZE]
x point of the track
TH2F * h_xvyvedepT[3]
Energy vs x vs y.
float z[maxSIZE]
z point of the track
TH2F * h_edepThres3
Energy deposited above threshold per time bin.
virtual void defineHisto() override
Defines the histograms.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:649
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.
STL namespace.