]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSPulseGenerator.cxx
Extract ALTRO sample generation to a separate class AliPHOSPulseGenerator
[u/mrichter/AliRoot.git] / PHOS / AliPHOSPulseGenerator.cxx
CommitLineData
431a9211 1/**************************************************************************
2 * Copyright(c) 2007, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16/* $Id$ */
17
18// The class which simulates the pulse shape from the PHOS FEE shaper,
19// make sampled amplitudes, digitize them.
20// Use case:
21// AliPHOSPulseGenerator *pulse = new AliPHOSPulseGenerator(energy,time);
22// Int_t *adcHG = new Int_t[pulse->GetRawFormatTimeBins()];
23// Int_t *adcLG= new Int_t[pulse->GetRawFormatTimeBins()];
24// pulse->MakeSamples();
25// pulse->GetSamples(adcHG, adcHG) ;
26// pulse->Print();
27// pulse->Draw();
28//
29// Author: Yuri Kharlov, after Yves Schutz and Per Thomas Hille
30
31// --- ROOT system ---
32#include "TMath.h"
33#include "TF1.h"
34#include "TGraph.h"
35#include "TCanvas.h"
36
37// --- AliRoot header files ---
38#include "AliLog.h"
39#include "AliPHOSPulseGenerator.h"
40
41// --- Standard library ---
42#include <cmath>
43#include <iostream>
44
45using std::cout;
46using std::endl;
47
48ClassImp(AliPHOSPulseGenerator)
49
50Double_t AliPHOSPulseGenerator::fgCapa = 1.; // 1pF
51Int_t AliPHOSPulseGenerator::fgOrder = 2 ; // order of the Gamma function
52Double_t AliPHOSPulseGenerator::fgTimeMax = 2.56E-5 ; // each sample is over 100 ns fTimeMax/fTimeBins
53Double_t AliPHOSPulseGenerator::fgTimePeak = 4.1E-6 ; // 4 micro seconds
54Double_t AliPHOSPulseGenerator::fgTimeTrigger = 100E-9 ; // 100ns, just for a reference
55Double_t AliPHOSPulseGenerator::fgHighCharge = 8.2; // adjusted for a high gain range of 5.12 GeV (10 bits)
56Double_t AliPHOSPulseGenerator::fgHighGain = 6.64;
57Double_t AliPHOSPulseGenerator::fgHighLowGainFactor = 16.; // adjusted for a low gain range of 82 GeV (10 bits)
58
59//-----------------------------------------------------------------------------
60AliPHOSPulseGenerator::AliPHOSPulseGenerator(Double_t a, Double_t t0)
61 : TObject(), fAmplitude(a), fTZero(t0), fDataHG(0), fDataLG(0), fDigitize(kTRUE)
62{
63 // Contruct a pulsegenrator object and initializes all necessary parameters
64 // @param a digit amplitude in GeV
65 // @param t0 time delay in nanoseconds of signal relative the first sample.
66 // This value should be between 0 and Ts, where Ts is the sample interval
67
68 fDataHG = new Double_t[fkTimeBins];
69 fDataLG = new Double_t[fkTimeBins];
70}
71
72//-----------------------------------------------------------------------------
73AliPHOSPulseGenerator::AliPHOSPulseGenerator(const AliPHOSPulseGenerator & pulse)
74 : TObject(), fAmplitude(pulse.fAmplitude), fTZero(pulse.fTZero), fDataHG(0), fDataLG(0), fDigitize(kTRUE)
75{
76 fDataHG = new Double_t[pulse.fkTimeBins];
77 fDataLG = new Double_t[pulse.fkTimeBins];
78 for (Int_t i=0; i<pulse.fkTimeBins; i++) {
79 fDataHG[i] = pulse.fDataHG[i];
80 fDataLG[i] = pulse.fDataHG[i];
81 }
82}
83
84//-----------------------------------------------------------------------------
85AliPHOSPulseGenerator::~AliPHOSPulseGenerator()
86{
87 // Destructor: delete arrays of samples
88
89 delete [] fDataHG;
90 fDataHG=0;
91 delete [] fDataLG;
92 fDataLG=0;
93}
94
95//-----------------------------------------------------------------------------
96void AliPHOSPulseGenerator::AddBaseline(Double_t baselineLevel)
97{
98 // Adds a baseline offset to the signal
99 // @param baselineLevel The basline level to add
100 for (Int_t i=0; i<fkTimeBins; i++) {
101 fDataHG[i] += baselineLevel;
102 fDataLG[i] += baselineLevel;
103 }
104}
105
106//-----------------------------------------------------------------------------
107void AliPHOSPulseGenerator::AddNoise(Double_t * /* sigma */)
108{
109 // Adds Gaussian white noise to the sample array given by *dataPtr.
110 // @param sigma the noise amplitude in entities of ADC levels
111
112 AliError("not implemented yet");
113}
114
115
116//-----------------------------------------------------------------------------
117void AliPHOSPulseGenerator::AddNoise(Double_t * /* sigma */, Double_t /* cutoff */)
118{
119 //Adds correlated Gaussian noise with cutof frequency "cutoff"
120 // @param sigma noise amplitude in entities of ADC levels
121 // @param -30DB cutoff frequency of the noise in entities of sampling frequency
122
123 AliError("not implemented yet");
124}
125
126//-----------------------------------------------------------------------------
127void AliPHOSPulseGenerator::AddPretriggerSamples(Int_t nPresamples)
128{
129 // Adds pretrigger samples to the sample arrays and replace them
130 // with concatinated and truncated arrays
131
132 Double_t *tmpDataHG = new Double_t[fkTimeBins];
133 Double_t *tmpDataLG = new Double_t[fkTimeBins];
134 Int_t i;
135 for (i=0; i<fkTimeBins; i++) {
136 tmpDataHG[i] = fDataHG[i];
137 tmpDataLG[i] = fDataLG[i];
138 }
139 for (i=0; i<fkTimeBins; i++) {
140 if (i<nPresamples) {
141 tmpDataHG[i] = 0.;
142 tmpDataLG[i] = 0.;
143 }
144 else {
145 tmpDataHG[i] = fDataHG[i-nPresamples];
146 tmpDataLG[i] = fDataLG[i-nPresamples];
147 }
148 }
149 delete [] tmpDataHG;
150 delete [] tmpDataLG;
151}
152
153//-----------------------------------------------------------------------------
154void AliPHOSPulseGenerator::Digitize()
155{
156 // Emulates ADC: rounds down to nearest integer value all amplitudes
157 for (Int_t i=0; i<fkTimeBins; i++) {
158 fDataHG[i] = (Double_t) ((Int_t)(fDataHG[i] + 0.5));
159 fDataLG[i] = (Double_t) ((Int_t)(fDataLG[i] + 0.5));
160 }
161}
162
163//-----------------------------------------------------------------------------
164Double_t AliPHOSPulseGenerator::RawResponseFunction(Double_t *x, Double_t *par)
165{
166 // Shape of the electronics raw reponse:
167 // It is a semi-gaussian, 2nd order Gamma function of the general form
168 // v(t) = n**n * Q * A**n / C *(t/tp)**n * exp(-n * t/tp) with
169 // tp : peaking time par[0]
170 // n : order of the function
171 // C : integrating capacitor in the preamplifier
172 // A : open loop gain of the preamplifier
173 // Q : the total APD charge to be measured Q = C * energy
174
175 Double_t signal ;
176 Double_t xx = x[0] - ( fgTimeTrigger + par[3] ) ;
177
178 if (xx < 0 || xx > fgTimeMax)
179 signal = 0. ;
180 else {
181 Double_t fac = par[0] * TMath::Power(fgOrder, fgOrder) * TMath::Power(par[1], fgOrder)/fgCapa ;
182 signal = fac * par[2] * TMath::Power(xx/fgTimePeak, fgOrder) * TMath::Exp(-fgOrder*(xx/fgTimePeak)) ;
183 }
184 return signal ;
185}
186
187//__________________________________________________________________
188Double_t AliPHOSPulseGenerator::RawResponseFunctionMax(Double_t charge, Double_t gain)
189{
190 // Maximum value of the shaper response function
191 return ( charge * TMath::Power(fgOrder, fgOrder) * TMath::Power(gain, fgOrder)
192 / ( fgCapa * TMath::Exp(fgOrder) ) );
193}
194
195//__________________________________________________________________
196Bool_t AliPHOSPulseGenerator::MakeSamples()
197{
198 // for a start time fTZero and an amplitude fAmplitude given by digit,
199 // calculates the raw sampled response AliPHOSPulseGenerator::RawResponseFunction
200
201 const Int_t kRawSignalOverflow = 0x3FF ;
202 Bool_t lowGain = kFALSE ;
203
204 TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
205
206 for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
207 signalF.SetParameter(0, fgHighCharge) ;
208 signalF.SetParameter(1, fgHighGain) ;
209 signalF.SetParameter(2, fAmplitude) ;
210 signalF.SetParameter(3, fTZero) ;
211 Double_t time = iTime * GetRawFormatTimeMax() / GetRawFormatTimeBins() ;
212 Double_t signal = signalF.Eval(time) ;
213 if ( static_cast<Int_t>(signal+0.5) > kRawSignalOverflow ){ // larger than 10 bits
214 signal = kRawSignalOverflow ;
215 lowGain = kTRUE ;
216 }
217 fDataHG[iTime] = signal;
218
219 signalF.SetParameter(0, GetRawFormatLowCharge() ) ;
220 signalF.SetParameter(1, GetRawFormatLowGain() ) ;
221 signal = signalF.Eval(time) ;
222 if ( static_cast<Int_t>(signal+0.5) > kRawSignalOverflow) // larger than 10 bits
223 signal = kRawSignalOverflow ;
224 fDataLG[iTime] = signal;
225
226 }
227 return lowGain ;
228}
229
230//__________________________________________________________________
231void AliPHOSPulseGenerator::GetSamples(Int_t *adcHG, Int_t *adcLG) const
232{
233 // Return integer sample arrays adcHG and adcLG
234 for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
235 adcHG[iTime] = static_cast<Int_t>(fDataHG[iTime]) ;
236 adcLG[iTime] = static_cast<Int_t>(fDataLG[iTime]) ;
237 }
238}
239
240//__________________________________________________________________
241void AliPHOSPulseGenerator::Print(Option_t*) const
242{
243 // Prints sampled amplitudes to stdout
244 Int_t i;
245 cout << "High gain: ";
246 for (i=0; i<fkTimeBins; i++)
247 cout << (Int_t)fDataHG[i] << " ";
248 cout << endl;
249
250 cout << "Low gain: ";
251 for (i=0; i<fkTimeBins; i++)
252 cout << (Int_t)fDataLG[i] << " ";
253 cout << endl;
254}
255
256//__________________________________________________________________
257void AliPHOSPulseGenerator::Draw(Option_t*)
258{
259 // Draw graphs with high and low gain samples
260
261 Double_t *time = new Double_t[fkTimeBins];
262 for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
263 time[iTime] = iTime * GetRawFormatTimeMax() / GetRawFormatTimeBins() ;
264 }
265 Int_t nPoints = fkTimeBins;
266 TGraph *graphHG = new TGraph(nPoints,time,fDataHG);
267 TGraph *graphLG = new TGraph(nPoints,time,fDataLG);
268 graphHG->SetMarkerStyle(20);
269 graphLG->SetMarkerStyle(24);
270 graphHG->SetMarkerSize(0.3);
271 graphLG->SetMarkerSize(0.3);
272 graphHG->SetTitle("High gain samples");
273 graphLG->SetTitle("Low gain samples");
274
275 TCanvas *c1 = new TCanvas("c1","Raw ALTRO samples",10,10,700,500);
276 c1->Divide(2,1);
277 c1->cd(1);
278 gPad->SetLeftMargin(0.15);
279 graphHG->Draw("AP");
280 graphHG->GetHistogram()->SetTitleOffset(1.6,"Y");
281 graphHG->GetHistogram()->SetXTitle("time, #musec");
282 graphHG->GetHistogram()->SetYTitle("Amplitude, ADC counts");
283 c1->cd(2);
284 gPad->SetLeftMargin(0.15);
285 graphLG->Draw("AP");
286 graphLG->GetHistogram()->SetTitleOffset(1.6,"Y");
287 graphLG->GetHistogram()->SetXTitle("time, #musec");
288 graphLG->GetHistogram()->SetYTitle("Amplitude, ADC counts");
289 c1->Update();
290}