Belle II Software  release-08-01-10
CsIStudyModule.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 <framework/core/HistoModule.h>
10 #include <framework/datastore/StoreArray.h>
11 #include <framework/datastore/StoreObjPtr.h>
12 #include <framework/dataobjects/EventMetaData.h>
13 #include <framework/gearbox/Const.h>
14 #include <beast/csi/modules/CsiModule.h>
15 #include <beast/csi/dataobjects/CsiSimHit.h>
16 #include <beast/csi/dataobjects/CsiHit.h>
17 #include <beast/csi/dataobjects/CsiDigiHit.h>
18 #include <beast/csi/modules/CsIStudyModule.h>
19 
20 // ROOT
21 #include <TH1.h>
22 #include <TFile.h>
23 
24 using namespace std;
25 using namespace Belle2;
26 using namespace csi;
27 
28 //-----------------------------------------------------------------
29 // Register the Module
30 //-----------------------------------------------------------------
31 REG_MODULE(CsIStudy);
32 
33 //-----------------------------------------------------------------
34 // Implementation
35 //-----------------------------------------------------------------
36 
37 CsIStudyModule::CsIStudyModule() : HistoModule(),
38  fWF(NULL),
39  h_CrystalEdep(NULL),
40  h_CrystalSpectrum(NULL),
41  h_CrystalRadDose(NULL),
42  h_CrystalRadDoseSH(NULL),
43  h_NhitCrystal(NULL),
44  h_LightYieldCrystal(NULL),
45  h_Waveform(NULL),
46  h_Gate(NULL),
47  h_Charge(NULL),
48  h_TrueEdep(NULL),
49  h_Height(NULL),
50  h_Weight(NULL)
51 
52 {
53  // Set module properties
54  setDescription("Analyze simulations of CsI readings in BEAST. Requires HistoManager module.");
55 
56  // Parameter definitions
57  addParam("nWFSamples", m_nWFSamples, "Number of samples in the saved waveforms", 8000);
58  addParam("paramTemplate", m_paramTemplate, "Template of an input parameter. Noop for now.", 0.0);
59 
60  string rootfilepath = "output/waveforms.root";
61  addParam("waveformFilename", m_waveformFilename, "Path to the root file to save waveforms", rootfilepath);
62 
63 }
64 
66 {
67 }
68 
70 {
71  h_CrystalEdep = new TH1F("Crystal_Edep", "Energy distribution in each crystal;CellID", 16, -0.5, 15.5);
72  h_CrystalSpectrum = new TH1F("Crystal_Edist", "Photon energy distribution in all crystals;", 100, 0, 0.1);
73  h_NhitCrystal = new TH1F("Crystal_HitRate", "Number of hit per crystal;CellID; hit/s", 16, -0.5, 15.5);
74  h_CrystalRadDose = new TH1F("Crystal_RadDose", "Crystal Radiation Dose;CellID;Gy/yr", 16, -0.5, 15.5);
75  h_CrystalRadDoseSH = new TH1F("Crystal_RadDose_SH", "Crystal Radiation Dose from SimHits;CellID;Gy/yr", 16, -0.5, 15.5);
76  h_LightYieldCrystal = new TH1F("Crystal_g_yield", "Light yield each crystal;CellID;gamma/sample", 16, -0.5, 15.5);
77  h_Waveform = new TH1S("Waveform", "Recorded waveform;Time index;ADC bits", m_nWFSamples, 0, m_nWFSamples - 1);
78  h_Gate = new TH1C("Gate", "Recorded gate;Time index;ADC bits", m_nWFSamples, 0, m_nWFSamples - 1);
79  h_Charge = new TH1F("Charge", "log_{10} of integrated Charge", 200, 0, 10);
80  h_TrueEdep = new TH1F("TrueEdep", "log_{10} of true deposited energy/GeV", 200, -6, 0);
81  h_Height = new TH1F("Height", "log10 of pulse height", 200, 0, 6);
82  h_Weight = new TH1F("Weight", "Distribution of the event weights form the generator", 1000, 0, 10);
83 }
84 
86 {
87  REG_HISTOGRAM // required to register histograms to HistoManager
88 
89  m_hits.isOptional();
90  m_simhits.isRequired();
91  m_digihits.isOptional();
92 
93  fWF = new TFile(m_waveformFilename.c_str(), "recreate");
94 }
95 
97 {
98 }
99 
101 {
102  // h_LightYieldCrystal->Divide( h_CrystalEdep );
103  // h_LightYieldCrystal->Scale( 1e-6 );
104 }
105 
106 
108 {
109 
110  StoreObjPtr<EventMetaData> eventMetaDataPtr;
111  double eventWeight = eventMetaDataPtr->getGeneratedWeight();
112 
113  h_Weight->Fill(eventWeight);
114 
115  //Loop over CsiHits
116  if (m_hits.getEntries() > 0) {
117  int hitNum = m_hits.getEntries();
119  double Mass = 5;
120  double edeptodose = GeVtoJ / Mass * usInYr / Sampletime;
123  for (int i = 0; i < hitNum; i++) {
124  CsiHit* aCsIHit = m_hits[i];
125  int m_cellID = aCsIHit->getCellId();
126  double edep = aCsIHit->getEnergyDep();
127  //double hitTime = aCsIHit->getTimeAve(); /**< Time of the hit*/
128 
129  // Fill histograms
130  h_CrystalSpectrum->Fill(edep, eventWeight);
131  h_CrystalEdep->Fill(m_cellID, edep * eventWeight);
132  h_CrystalRadDose->Fill(m_cellID, edep * edeptodose * eventWeight);
133 
134  //Number of hits per second
135  h_NhitCrystal->Fill(m_cellID, eventWeight * 1.0e9 / Sampletime);
136  }
137  }
138 
139  //Loop over CsiSimHits
140  if (m_simhits.getEntries() > 0) {
141  int hitNum = m_simhits.getEntries();
143  double Mass = 5;
144  double edeptodose = GeVtoJ / Mass * usInYr / Sampletime;
147  for (int i = 0; i < hitNum; i++) {
148  CsiSimHit* aCsIHit = m_simhits[i];
149  int m_cellID = aCsIHit->getCellId();
150  double edep = aCsIHit->getEnergyDep();
152  // Fill histograms
153  h_CrystalRadDoseSH->Fill(m_cellID, eventWeight * edep * edeptodose);
154 
155  // To get the Number of photons per GeV (divide by total edep before plotting)
156  if (Const::photon.getPDGCode() == aCsIHit->getPDGCode()) {
157  h_LightYieldCrystal->Fill(m_cellID, eventWeight);
158  }
159  }
160  }
161 
162  //Loop over Digihits
163  if (m_digihits.getEntries() > 0) {
164  int hitNum = m_digihits.getEntries();
167  for (int i = 0; i < hitNum; i++) {
168  CsiDigiHit* aCsIDigiHit = m_digihits[i];
169 
170  int cellID = aCsIDigiHit->getCellId();
171  // cast to double before taking the log
172  double charge = (double) aCsIDigiHit->getCharge();
173  double height = (double) aCsIDigiHit->getMaxVal();
174  double baseline = (double) aCsIDigiHit->getBaseline();
175  double trueEdep = aCsIDigiHit->getTrueEdep();
176  vector<uint16_t>* waveform = aCsIDigiHit->getWaveform();
177  vector<uint8_t>* status = aCsIDigiHit->getStatusBits();
178 
179  h_TrueEdep->Fill(log10(trueEdep), eventWeight);
180  h_Charge->Fill(log10(charge - baseline), eventWeight);
181  h_Height->Fill(log10(height), eventWeight);
182 
183  // Write the wavforms. All have the same number in the TFile, but
184  // we can access them by their cycle number:
185  // root[] TFile *f1 = new TFile("filename.root")
186  // root[] TH1F *Waveform4; f1->GetObject("h_Waveform;4",Waveform4);
187 
188  if (waveform->size()) {
189  h_Gate->Reset();
190  h_Waveform->Reset();
191 
192  char histoTitle[80];
193  sprintf(histoTitle, "Waveform Hit No. %i, Cell No %i", i, cellID);
194  h_Waveform->SetTitle(histoTitle);
195 
196 
197  for (uint iBin = 0; iBin < waveform->size(); iBin++) {
198  bool* statusBits = readDPPStatusBits(status->at(iBin));
199  h_Gate->Fill(iBin + 1, (int) statusBits[1]);
200  h_Waveform->Fill(iBin + 1, waveform->at(iBin));
201  }
202 
203  fWF->cd();
204  h_Gate->Write();
205  h_Waveform->Write();
206  }
207  }
208  }
209 }
210 
211 
213 {
214  fWF->Close();
215 }
216 
217 
219 {
220 
221  bool* bit = new bool[8];
222 
223  // fill data
224  for (int i = 0; i < 8; i++) {
225  bit[i] = ((data >> i) & 0x01);
226  }
227 
228  return bit;
229 }
static const ParticleType photon
photon particle
Definition: Const.h:664
const double Sampletime
Sample time in us.
TH1F * h_CrystalEdep
Energy deposited in each crystal.
TH1F * h_Charge
Distribution of the integrated charge from the pulse processing algorithm.
TH1F * h_CrystalSpectrum
Distribution of the energy of each simhit.
TH1F * h_TrueEdep
Distribution of the true total deposited energy in the event.
virtual ~CsIStudyModule()
Empty destructor.
virtual void initialize() override
Initialization: building histograms.
virtual void event() override
Read each event, calculate doses and fill histograms.
bool * readDPPStatusBits(char packedbits)
Reads and unpacks the status bits of the digitizer Unpacks in the appropriate booleans.
StoreArray< CsiHit > m_hits
Array of the crystal hits (1 per event)
double m_paramTemplate
Template of an input parameter.
virtual void endRun() override
To do at the end of each runs.
TH1F * h_Height
Max height of the peak signa(.
virtual void terminate() override
Clean everything up.
StoreArray< CsiSimHit > m_simhits
Array of the simulated hits in the crystals (incl.
int m_nWFSamples
Number of samples in each of the waveforms (in the near future this should be imported from the digi ...
TH1F * h_Weight
Event weight assigned by the generator.
TH1S * h_Waveform
Contains the digitized waveform (if enabled in the digitizer)
const double usInYr
Conversion between us to nominal year.
virtual void beginRun() override
To do at the beginning of each runs.
TH1F * h_LightYieldCrystal
Number of photons hits in each crystal (to validate light yield..)
TH1C * h_Gate
Contains the integration gate (if waveforms are enabled in the digitizer)
TH1F * h_NhitCrystal
Number of hits in each crystal (to get hit rate)
StoreArray< CsiDigiHit > m_digihits
Array of the crystal digitized hits (1 per event per fired crystal)
TH1F * h_CrystalRadDoseSH
Yearly radiation dose deposited in each crystal (for hit-simhit check)
std::string m_waveformFilename
Path where to save the waveforms (root file)
const double GeVtoJ
Conversion between GeV to Joule.
TFile * fWF
ROOT file to save waveforms.
TH1F * h_CrystalRadDose
Yearly radiation dose deposited in each crystal.
virtual void defineHisto() override
function to define histograms
Class to store Csi digitized hits (output of CsIDigitizer) relation to CsiHit filled in beast/csi/mod...
Definition: CsiDigiHit.h:26
double getTrueEdep() const
Get true deposited energy.
Definition: CsiDigiHit.h:77
std::vector< uint8_t > * getStatusBits()
Get Status bits array.
Definition: CsiDigiHit.h:114
uint16_t getBaseline() const
Get Baseline.
Definition: CsiDigiHit.h:87
uint32_t getCharge() const
Get Integrated Charge.
Definition: CsiDigiHit.h:67
std::vector< uint16_t > * getWaveform()
Get Waveform array.
Definition: CsiDigiHit.h:105
uint8_t getCellId() const
Get Cell ID.
Definition: CsiDigiHit.h:61
uint16_t getMaxVal() const
Get Maximal Value.
Definition: CsiDigiHit.h:95
Class to store simulated hits which equate to average of CsiSimHit on crystals Input for digitization...
Definition: CsiHit.h:26
int getCellId() const
Get Cell ID.
Definition: CsiHit.h:73
double getEnergyDep() const
Get deposit energy.
Definition: CsiHit.h:78
ClassCsiSimHit - Geant4 simulated hits in CsI crystals in BEAST.
Definition: CsiSimHit.h:31
int getPDGCode() const
Get Particle PDG (can be one of secondaries)
Definition: CsiSimHit.h:97
int getCellId() const
Get Cell ID.
Definition: CsiSimHit.h:87
double getEnergyDep() const
Get Deposit energy.
Definition: CsiSimHit.h:107
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
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.