Bartender LYSO
Digitizer simulation for the LYSO calorimeter prototype
Loading...
Searching...
No Matches
bar.cc
Go to the documentation of this file.
1
5#include "bar.hh"
6
7using namespace std;
8using namespace TMath;
9
10
11BarLYSO::BarLYSO(const char*inputFilename, Int_t threadID)
12{
13 // Set the threadID
14 fThreadID = threadID;
15
16 // Initialize the random generators
17 fRandPars = new TRandom3(0);
18 fRandNoise = new TRandom3(0);
19
20 // Initialize DAQ
21 fDAQ = new DAQ();
22
23 // Initialize whole WF containers [CHANNELS]x[SAMPLINGS]
24 fEvent = -1;
25 fFront.resize(CHANNELS, vector<Float_t>(SAMPLINGS, 0));
26 fBack.resize(CHANNELS, vector<Float_t>(SAMPLINGS, 0));
27 fTimes_F.resize(CHANNELS, vector<Float_t>(SAMPLINGS, 0));
28 fTimes_B.resize(CHANNELS, vector<Float_t>(SAMPLINGS, 0));
29
30 // Determine the output filename based on the BarLYSO ID
31 string outputFilename = GenerateOutputFilename(inputFilename);
32
33 // Create the file.root and the TTrees
34 fOutFile = TFile::Open(outputFilename.c_str(), "RECREATE");
35
36 fOutTree = new TTree("lyso_wfs", "lyso_wfs");
37 fOutTree->Branch("Event", &fEvent);
38 fOutTree->Branch("Front", &fFront);
39 fOutTree->Branch("Back", &fBack);
40
41 fTimesTree = new TTree("lyso_wfs_times", "lyso_wfs_times");
42 fTimesTree->Branch("Time_F", &fTimes_F);
43 fTimesTree->Branch("Time_B", &fTimes_B);
44}
45
46
47
49{
50 // Delete DAQ and random generators
51 delete fDAQ;
52 delete hPars;
53 delete fRandPars;
54 delete fRandNoise;
55
56 // Ensure that we delete the TTree and TFile objects only if they are not null
57 if(fOutTree)
58 {
59 fOutTree->ResetBranchAddresses(); // Reset any attached branches if needed
60 delete fOutTree;
61 }
62 if(fTimesTree)
63 {
64 fTimesTree->ResetBranchAddresses(); // Reset any attached branches if needed
65 delete fTimesTree;
66 }
67 if(fOutFile)
68 {
69 fOutFile->Close(); // Make sure to close the file properly
70 delete fOutFile;
71 }
72}
73
74
75
76string BarLYSO::GenerateOutputFilename(const char* inputFilename)
77{
78 // Regular expression to extract digits following "MCID_"
79 regex regex("\\bMCID_(\\d+)");
80
81 // Storage for matched results
82 smatch match;
83
84 // Convert inputFilename to string
85 string filename = inputFilename;
86
87 // Check if there is a match for digits following "MCID_" in the filename
88 if(regex_search(filename, match, regex))
89 {
90 // Extract the matched digits as a string
91 string runID = match[1].str();
92
93 if (!runID.empty())
94 {
95 // Convert the extracted string to an integer and assign it to the fID member variable
96 fID = stoi(runID);
97 }
98 }
99
100 // Determine the output filename based on the Bar ID
101 string outputFilename;
102 if(fID)
103 {
104 // If fThreadID is -1, omit the thread ID part in the filename
105 if(fThreadID == -1)
106 {
107 outputFilename = "./RootFiles/BarID_" + to_string(fID) + ".root";
108 }
109 else
110 {
111 outputFilename = "./RootFiles/BarID_" + to_string(fID) + "_t" + to_string(fThreadID) + ".root";
112 }
113 }
114 else
115 {
116 outputFilename = "./RootFiles/output.root";
117 }
118
119 return outputFilename;
120}
121
122
123
124void BarLYSO::SetSamplingTimes()
125{
126 // Set the times
127 for(Int_t j = 0; j < CHANNELS; j++)
128 {
129 if(fDAQ->fIsBinSizeConstant)
130 {
131 for(Int_t i = 0; i < SAMPLINGS; i++)
132 {
133 fTimes_F[j][i] = (Float_t) i / fDAQ->fSamplingSpeed;
134 fTimes_B[j][i] = (Float_t) i / fDAQ->fSamplingSpeed;
135 }
136 }
137 else
138 {
139 // Initialize the first time points
140 fTimes_F[j][0] = 0.0;
141 fTimes_B[j][0] = 0.0;
142
143 for(Int_t i = 0; i < SAMPLINGS - 1; i++) // Up to SAMPLINGS - 1
144 {
145 Float_t bin_F, bin_B;
146 do
147 {
148 bin_F = fDAQ->binRand->Gaus(1.0 / fDAQ->fSamplingSpeed, fDAQ->fSigmaBinSize);
149 } while (bin_F < 0.5 * (1.0 / fDAQ->fSamplingSpeed) || bin_F > 1.5 * (1.0 / fDAQ->fSamplingSpeed));
150
151 do
152 {
153 bin_B = fDAQ->binRand->Gaus(1.0 / fDAQ->fSamplingSpeed, fDAQ->fSigmaBinSize);
154 } while (bin_B < 0.5 * (1.0 / fDAQ->fSamplingSpeed) || bin_B > 1.5 * (1.0 / fDAQ->fSamplingSpeed));
155
156 fTimes_F[j][i+1] = fTimes_F[j][i] + bin_F;
157 fTimes_B[j][i+1] = fTimes_B[j][i] + bin_B;
158 }
159 }
160 }
161
162 fTimesTree->Fill();
163}
164
165
166
168{
169 Int_t status;
170 Double_t charge, A, tau_rise, tau_dec;
171
172 // Get a TTree from the text file of parameters
173 TTree *tree = new TTree("tree", "mytree");
174 tree->ReadFile(fInputFilename.c_str());
175
176 tree->SetBranchAddress("Status", &status);
177 tree->SetBranchAddress("Charge", &charge);
178 tree->SetBranchAddress("A", &A);
179 tree->SetBranchAddress("Tau_rise", &tau_rise);
180 tree->SetBranchAddress("Tau_fall", &tau_dec);
181
182 // Create hPars
183 hPars = new TH3D("hPars", "Fitted All", fHisto_A[0], fHisto_A[1], fHisto_A[2], fHisto_Tau_rise[0], fHisto_Tau_rise[1], fHisto_Tau_rise[2], fHisto_Tau_dec[0], fHisto_Tau_dec[1], fHisto_Tau_dec[2]);
184
185 // Fill hPars
186 for(Int_t i = 0; i < tree->GetEntries(); i++)
187 {
188 tree->GetEntry(i);
189
190 // Check conditions about charge
191 if(status==0 && charge >= fChargeCuts[0] && charge <= fChargeCuts[1])
192 {
193 hPars->Fill(A, tau_rise, tau_dec);
194 }
195 }
196
197 // Delete the TTree
198 delete tree;
199}
200
201
202
203Float_t BarLYSO::Wave_OnePhel(Float_t t, Double_t A, Double_t tau_rise, Double_t tau_dec, Double_t timePhel)
204{
205 Float_t funcVal;
206 Double_t expRise = Exp(-(t-timePhel)/tau_rise);
207 Double_t expDec = Exp(-(t-timePhel)/tau_dec);
208
209 // 1-Phel function
210 funcVal = static_cast<Float_t>(A*(expRise-expDec)*((t > timePhel) ? 1:0));
211
212 // Numerical fixing
213 if(funcVal > 0 || IsNaN(funcVal)) funcVal = 0;
214
215 return funcVal;
216}
217
218
219
221{
222 fEvent = event;
223 fFront.assign(CHANNELS, std::vector<Float_t>(SAMPLINGS, 0.));
224 fBack.assign(CHANNELS, std::vector<Float_t>(SAMPLINGS, 0.));
225}
226
227
228
230{
231 for(Int_t ch = 0; ch < CHANNELS; ch++)
232 {
233 fFront[ch].clear();
234 fBack[ch].clear();
235 }
236
237 fFront.clear();
238 fBack.clear();
239}
240
241
242
243void BarLYSO::SetFrontWaveform(Int_t channel, Double_t start)
244{
245 // Sample from hPars the parameters of 1-Phel WF
246 Double_t A, tau_rise, tau_dec;
247 hPars->GetRandom3(A, tau_rise, tau_dec, fRandPars);
248
249 // Evaluate and sum the new 1-Phel WF to the existing one
250 for(Int_t bin = 0; bin < SAMPLINGS; bin++)
251 {
252 fFront[channel][bin] += Wave_OnePhel(fTimes_F[channel][bin], A, tau_rise, tau_dec, start + ZERO_TIME_BIN);
253 }
254}
255
256
257
258void BarLYSO::SetBackWaveform(Int_t channel, Double_t start)
259{
260 // Sample from hPars the parameters of 1-Phel WF
261 Double_t A, tau_rise, tau_dec;
262 hPars->GetRandom3(A, tau_rise, tau_dec, fRandPars);
263
264 // Evaluate and sum the new 1-Phel WF to the existing one
265 for(Int_t bin = 0; bin < SAMPLINGS; bin++)
266 {
267 fBack[channel][bin] += Wave_OnePhel(fTimes_B[channel][bin], A, tau_rise, tau_dec, start + ZERO_TIME_BIN);
268 }
269}
270
271
272
274{
275 // Recompute the gain and add noise
276 Float_t k = fDAQ->ComputeFactorOfGainConversion();
277 for(Int_t ch = 0; ch < CHANNELS; ch++)
278 {
279 for(Int_t bin = 0; bin < SAMPLINGS; bin++)
280 {
281 // Set every bin with gaussian noise
282 fFront[ch][bin] = k*fFront[ch][bin] + Add_Noise();
283 fBack[ch][bin] = k*fBack[ch][bin] + Add_Noise();
284 }
285 }
286
287 fOutTree->Fill();
288}
289
290
291
292void BarLYSO::SaveBar()
293{
294 fOutFile->cd();
295 fOutTree->Write("lyso_wfs");
296 fTimesTree->Write("lyso_wfs_times");
297}
Declaration of the class BarLYSO.
std::vector< std::vector< Float_t > > fBack
Container for Back-Detector waveforms: a 3-dimensional matrix with indices for event,...
Definition bar.hh:165
Float_t Wave_OnePhel(Float_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:203
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
Double_t fChargeCuts[2]
Cuts in the charge spectrum of input best fit parameters data; [0] represents the minimum,...
Definition bar.hh:177
Float_t Add_Noise()
Returns the value of the noise.
Definition bar.hh:198
void SaveEvent()
Saves all the samples from fFront and fBack into a text file.
Definition bar.cc:273
void InitializeBaselines(Int_t event)
Method to initialize the entire fFront and fBack with a noise baseline.
Definition bar.cc:220
void SetParsDistro()
Sets the 3D histogram hPars.
Definition bar.cc:167
std::vector< std::vector< Float_t > > fFront
Container for Front-Detector waveforms: a 3-dimensional matrix with indices for event,...
Definition bar.hh:164
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:258
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 ClearContainers()
Method to add a I-Phel waveform to the corresponding event and channel of the Front-Detector.
Definition bar.cc:229
~BarLYSO()
Destructor of the class.
Definition bar.cc:48
TRandom3 * fRandPars
Random generator for SetFrontWaveform() and SetBackWaveform()
Definition bar.hh:171
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 fID
Run ID of the Monte Carlo.
Definition bar.hh:160
Definition daq.hh:12
constexpr Int_t SAMPLINGS
Number of samplings for one waveform.
Definition globals.hh:10
constexpr Float_t ZERO_TIME_BIN
Delay of all waveforms in the [0, 1023] bins window.
Definition globals.hh:12
constexpr Int_t CHANNELS
Number of channels of the detectors.
Definition globals.hh:9