Bartender LYSO
Digitizer simulation for the LYSO calorimeter prototype
Loading...
Searching...
No Matches
bar.hh
Go to the documentation of this file.
1
5#ifndef BAR_HH
6#define BAR_HH
7
8#include <iostream>
9#include <fstream>
10#include <vector>
11#include <cfloat>
12#include <string>
13#include <regex>
14
15#include <TH3D.h>
16#include <TRandom3.h>
17#include <TTree.h>
18#include <TMath.h>
19#include <TFile.h>
20
21#include "globals.hh"
22#include "daq.hh"
23
28{
29public:
37 BarLYSO(const char* inputFilename, Int_t threadID);
41 ~BarLYSO();
42
54 void SetParsDistro();
65 void InitializeBaselines(Int_t event);
81 void ClearContainers();
82 void SetSamplingTimes();
83 void SetFrontWaveform(Int_t channel, Double_t start);
97 void SetBackWaveform(Int_t channel, Double_t start);
108 void SaveEvent();
109 void SaveBar();
110
111 inline void SetSigmaNoise(Double_t newSigmaNoise) { fSigmaNoise = newSigmaNoise; }
112 inline void SetInputFilename(std::string newInputFilename) { fInputFilename = newInputFilename; }
117 inline void SetChargeCuts(Double_t min, Double_t max)
118 {
119 fChargeCuts[0] = min;
120 fChargeCuts[1] = max;
121 }
126 inline void SetHisto_A(Double_t nbins, Double_t min, Double_t max)
127 {
128 fHisto_A[0] = TMath::Nint(nbins);
129 fHisto_A[1] = min;
130 fHisto_A[2] = max;
131 }
136 inline void SetHisto_Tau_rise(Double_t nbins, Double_t min, Double_t max)
137 {
138 fHisto_Tau_rise[0] = TMath::Nint(nbins);
139 fHisto_Tau_rise[1] = min;
140 fHisto_Tau_rise[2] = max;
141 }
146 inline void SetHisto_Tau_dec(Double_t nbins, Double_t min, Double_t max)
147 {
148 fHisto_Tau_dec[0] = TMath::Nint(nbins);
149 fHisto_Tau_dec[1] = min;
150 fHisto_Tau_dec[2] = max;
151 }
152
153 inline void SetEvents(Int_t events) { EVENTS = events; }
154 inline Int_t GetEvents() const { return EVENTS; }
155 inline Int_t GetID() const { return fID; }
156 inline DAQ *GetDAQ() const { return fDAQ; }
157
158private:
159 Int_t fEvent;
160 Int_t fID;
161 Int_t EVENTS;
162 Int_t fThreadID;
163
164 std::vector<std::vector<Double_t>> fFront;
165 std::vector<std::vector<Double_t>> fBack;
166 std::vector<std::vector<Double_t>> fTimes_F;
167 std::vector<std::vector<Double_t>> fTimes_B;
168
169 TH3D *hPars;
171 TRandom3 *fRandPars;
172 TRandom3 *fRandNoise;
174 Double_t fSigmaNoise;
176 std::string fInputFilename;
177 Double_t fChargeCuts[2];
178 Double_t fHisto_A[3];
179 Double_t fHisto_Tau_rise[3];
180 Double_t fHisto_Tau_dec[3];
182 DAQ *fDAQ;
183
184 std::string GenerateOutputFilename(const char *inputFilename);
185
186 TFile *fOutFile = nullptr;
187 TTree *fOutTree = nullptr;
188
197 inline Double_t Add_Noise() { return fRandNoise->Gaus(BASELINE, fDAQ->fSigmaNoise); };
205 Double_t Wave_OnePhel(Double_t t, Double_t A, Double_t tau_rise, Double_t tau_dec, Double_t timePhel);
206
207};
208
209#endif // BAR_HH
Class for managing waveform construction for all events and channels.
Definition bar.hh:28
Double_t fHisto_Tau_dec[3]
Settings for histogram hPars related to parameter Tau_dec: [0] for number of bins,...
Definition bar.hh:180
TRandom3 * fRandNoise
Random generator for Add_Noise()
Definition bar.hh:172
Int_t fEvent
Number of events in the run.
Definition bar.hh:159
std::string fInputFilename
Name of the txt file of the best fit parameters data. See the introduction for more details about the...
Definition bar.hh:176
std::vector< std::vector< Double_t > > fBack
Container for Back-Detector waveforms: a 3-dimensional matrix with indices for event,...
Definition bar.hh:165
Double_t fChargeCuts[2]
Cuts in the charge spectrum of input best fit parameters data; [0] represents the minimum,...
Definition bar.hh:177
void SetSigmaNoise(Double_t newSigmaNoise)
Set fSigmaNoise, the noise of the DAQ.
Definition bar.hh:111
void SaveEvent()
Saves all the samples from fFront and fBack into a text file.
Definition bar.cc:263
void InitializeBaselines(Int_t event)
Method to initialize the entire fFront and fBack with a noise baseline.
Definition bar.cc:210
void SetInputFilename(std::string newInputFilename)
Set the name of the text file of the best fit parameters data.
Definition bar.hh:112
void SetParsDistro()
Sets the 3D histogram hPars.
Definition bar.cc:157
Double_t Add_Noise()
Returns the value of the noise.
Definition bar.hh:197
BarLYSO(const char *inputFilename, Int_t threadID)
Constructor of the class.
Definition bar.cc:11
void SetBackWaveform(Int_t channel, Double_t start)
Method to add a I-Phel waveform to the corresponding event and channel of the Back-Detector.
Definition bar.cc:248
TH3D * hPars
3D Histogram of One-Phel waveform parameters from which sampling will occur
Definition bar.hh:169
Double_t fHisto_A[3]
Settings for histogram hPars related to parameter A: [0] for number of bins, [1] for the lower limit,...
Definition bar.hh:178
void SetHisto_Tau_rise(Double_t nbins, Double_t min, Double_t max)
Set the number of bins, the lower and the upper limit related to parameter Tau_rise for the histogram...
Definition bar.hh:136
void ClearContainers()
Method to add a I-Phel waveform to the corresponding event and channel of the Front-Detector.
Definition bar.cc:219
~BarLYSO()
Destructor of the class.
Definition bar.cc:45
void SetHisto_A(Double_t nbins, Double_t min, Double_t max)
Set the number of bins, the lower and the upper limit related to parameter A for the histogram hPars.
Definition bar.hh:126
Double_t Wave_OnePhel(Double_t t, Double_t A, Double_t tau_rise, Double_t tau_dec, Double_t timePhel)
Returns the value at t of the analytical form of the One Photo-Electron waveform.
Definition bar.cc:193
TRandom3 * fRandPars
Random generator for SetFrontWaveform() and SetBackWaveform()
Definition bar.hh:171
Int_t GetID() const
Returns the ID of the Monte Carlo.
Definition bar.hh:155
Double_t fHisto_Tau_rise[3]
Settings for histogram hPars related to parameter Tau_rise: [0] for number of bins,...
Definition bar.hh:179
Int_t GetEvents() const
Returns the number of events in the run.
Definition bar.hh:154
Double_t fSigmaNoise
Noise of the DAQ, evaluated as the stDev of the pedestal distribution.
Definition bar.hh:174
void SetHisto_Tau_dec(Double_t nbins, Double_t min, Double_t max)
Set the number of bins, the lower and the upper limit related to parameter Tau_dec for the histogram ...
Definition bar.hh:146
void SetChargeCuts(Double_t min, Double_t max)
Set the cuts in the charge spectrum of input best fit parameters data.
Definition bar.hh:117
Int_t fID
Run ID of the Monte Carlo.
Definition bar.hh:160
std::vector< std::vector< Double_t > > fFront
Container for Front-Detector waveforms: a 3-dimensional matrix with indices for event,...
Definition bar.hh:164
Definition daq.hh:12
Double_t fSigmaNoise
Noise of the DAQ, evaluated as the stDev of the pedestal distribution.
Definition daq.hh:33
Definition of global constant values for the simulation.