]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSPulseGenerator.cxx
New class AliPHOSRawDecoder introduced.
[u/mrichter/AliRoot.git] / PHOS / AliPHOSPulseGenerator.cxx
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->AddNoise(1.);
25 //   pulse->MakeSamples();
26 //   pulse->GetSamples(adcHG, adcHG) ; 
27 //   pulse->Print();
28 //   pulse->Draw();
29 //
30 // Author: Yuri Kharlov, after Yves Schutz and Per Thomas Hille
31
32 // --- ROOT system ---
33 #include "TMath.h" 
34 #include "TF1.h" 
35 #include "TGraph.h" 
36 #include "TCanvas.h" 
37 #include "TRandom.h"
38
39 // --- AliRoot header files ---
40 #include "AliLog.h"
41 #include "AliPHOSPulseGenerator.h"
42
43 // --- Standard library ---
44 #include <cmath>
45 #include <iostream>
46
47 using std::cout;
48 using std::endl; 
49
50 ClassImp(AliPHOSPulseGenerator) 
51
52 Double_t AliPHOSPulseGenerator::fgCapa        = 1.1;       // 1pF 
53 Int_t    AliPHOSPulseGenerator::fgOrder       = 2 ;        // order of the Gamma function
54 Double_t AliPHOSPulseGenerator::fgTimePeak    = 2.1E-6 ;   // tau=2.1 micro seconds
55 Double_t AliPHOSPulseGenerator::fgTimeTrigger = 100E-9 ;   // one tick 100 ns
56 Double_t AliPHOSPulseGenerator::fgHighCharge  = 8.8;       // adjusted for a high gain range of 5 GeV (10 bits)
57 Double_t AliPHOSPulseGenerator::fgHighGain    = 6.85;
58 Double_t AliPHOSPulseGenerator::fgHighLowGainFactor = 16.; // adjusted for a low gain range of 80 GeV (10 bits) 
59
60 //-----------------------------------------------------------------------------
61 AliPHOSPulseGenerator::AliPHOSPulseGenerator(Double_t a, Double_t t0)
62   : TObject(), fAmplitude(a), fTZero(t0), fDataHG(0), fDataLG(0), fDigitize(kTRUE)
63 {
64   // Contruct a pulsegenrator object and initializes all necessary parameters
65   // @param a digit amplitude in GeV
66   // @param t0 time delay in nanoseconds of signal relative the first sample. 
67   // This value should be between 0 and Ts, where Ts is the sample interval
68
69   fDataHG = new Double_t[fkTimeBins];
70   fDataLG = new Double_t[fkTimeBins];
71   Reset();
72 }
73
74 //-----------------------------------------------------------------------------
75 AliPHOSPulseGenerator::AliPHOSPulseGenerator(const AliPHOSPulseGenerator & pulse)
76   : TObject(), fAmplitude(pulse.fAmplitude), fTZero(pulse.fTZero), fDataHG(0), fDataLG(0), fDigitize(kTRUE)
77 {
78   fDataHG = new Double_t[pulse.fkTimeBins];
79   fDataLG = new Double_t[pulse.fkTimeBins];
80   for (Int_t i=0; i<pulse.fkTimeBins; i++) {
81     fDataHG[i] = pulse.fDataHG[i];
82     fDataLG[i] = pulse.fDataHG[i];
83   }
84 }
85
86 //-----------------------------------------------------------------------------
87 AliPHOSPulseGenerator::~AliPHOSPulseGenerator()
88 {
89   // Destructor: delete arrays of samples
90
91   delete [] fDataHG;
92   fDataHG=0;
93   delete [] fDataLG;
94   fDataLG=0;
95 }
96
97 //-----------------------------------------------------------------------------
98 void AliPHOSPulseGenerator::Reset()
99 {
100   // Reset all sample amplitudes to 0
101
102   for (Int_t i=0; i<fkTimeBins; i++) {
103     fDataHG[i] = 0.;
104     fDataLG[i] = 0.;
105   }
106 }
107
108 //-----------------------------------------------------------------------------
109 void AliPHOSPulseGenerator::AddBaseline(Double_t baselineLevel)
110 {
111   // Adds a baseline offset to the signal
112   // @param baselineLevel The basline level to add
113   for (Int_t i=0; i<fkTimeBins; i++) {
114     fDataHG[i] += baselineLevel;
115     fDataLG[i] += baselineLevel;
116   }
117   // Digitize floating point amplitudes to integers
118   if (fDigitize) Digitize();
119 }
120
121 //-----------------------------------------------------------------------------
122 void AliPHOSPulseGenerator::AddNoise(Double_t sigma)
123 {
124   // Adds Gaussian uncorrelated to the sample array
125   // @param sigma the noise amplitude in entities of ADC levels  
126   
127   for (Int_t i=0; i<fkTimeBins; i++) {
128     fDataHG[i] = gRandom->Gaus(0., sigma) ; 
129     fDataLG[i] = gRandom->Gaus(0., sigma) ; 
130   }
131 }
132
133 //-----------------------------------------------------------------------------
134 void AliPHOSPulseGenerator::AddNoise(Double_t * /* sigma */, Double_t /* cutoff */)
135 {
136   //Adds correlated Gaussian noise with cutof frequency "cutoff"
137   // @param sigma noise amplitude in entities of ADC levels
138   // @param -30DB cutoff frequency of the noise in entities of sampling frequency
139
140   AliError("not implemented yet");
141 }
142
143 //-----------------------------------------------------------------------------
144 void AliPHOSPulseGenerator::AddPretriggerSamples(Int_t nPresamples)
145 {
146   // Adds pretrigger samples to the sample arrays and replace them
147   // with concatinated and truncated arrays
148
149   Double_t *tmpDataHG = new Double_t[fkTimeBins];
150   Double_t *tmpDataLG = new Double_t[fkTimeBins];
151   Int_t i;
152   for (i=0; i<fkTimeBins; i++) {
153     tmpDataHG[i] = fDataHG[i];
154     tmpDataLG[i] = fDataLG[i];
155   }
156   for (i=0; i<fkTimeBins; i++) {
157     if (i<nPresamples) {
158       fDataHG[i] = 0.;
159       fDataLG[i] = 0.;
160     }
161     else {
162       fDataHG[i] = tmpDataHG[i-nPresamples];
163       fDataLG[i] = tmpDataLG[i-nPresamples];
164     }
165   }
166   delete [] tmpDataHG;
167   delete [] tmpDataLG;
168 }
169
170 //-----------------------------------------------------------------------------
171 void AliPHOSPulseGenerator::Digitize()
172 {
173   // Emulates ADC: rounds down to nearest integer value all amplitudes
174   for (Int_t i=0; i<fkTimeBins; i++) {
175     fDataHG[i] = (Double_t) ((Int_t)(fDataHG[i] + 0.5));
176     fDataLG[i] = (Double_t) ((Int_t)(fDataLG[i] + 0.5));
177   }
178 }
179
180 //-----------------------------------------------------------------------------
181 Double_t AliPHOSPulseGenerator::RawResponseFunction(Double_t *x, Double_t *par) 
182 {
183   // Shape of the electronics raw reponse:
184   // It is a semi-gaussian, 2nd order Gamma function of the general form
185   // v(t) = n**n * Q * A**n / C *(t/tp)**n * exp(-n * t/tp) with 
186   // tp : peaking time par[0]
187   // n  : order of the function
188   // C  : integrating capacitor in the preamplifier
189   // A  : open loop gain of the preamplifier
190   // Q  : the total APD charge to be measured Q = C * energy
191   
192   Double_t signal ;
193   Double_t xx = x[0] - ( fgTimeTrigger + par[3] ) ; 
194
195   if (xx < 0 || xx > GetRawFormatTimeMax()) 
196     signal = 0. ;  
197   else {
198     Double_t fac = par[0] * TMath::Power(fgOrder, fgOrder) * TMath::Power(par[1], fgOrder)/fgCapa ; 
199     signal = fac * par[2] * TMath::Power(xx/fgTimePeak, fgOrder) * TMath::Exp(-fgOrder*(xx/fgTimePeak)) ;
200   }
201   return signal ;  
202 }
203
204 //__________________________________________________________________
205 Double_t AliPHOSPulseGenerator::RawResponseFunctionMax(Double_t charge, Double_t gain) 
206 {
207   // Maximum value of the shaper response function
208   return ( charge * TMath::Power(fgOrder, fgOrder) * TMath::Power(gain, fgOrder) 
209            / ( fgCapa * TMath::Exp(fgOrder) ) );  
210 }
211
212 //__________________________________________________________________
213 Bool_t AliPHOSPulseGenerator::MakeSamples()
214 {
215   // for a start time fTZero and an amplitude fAmplitude given by digit, 
216   // calculates the raw sampled response AliPHOSPulseGenerator::RawResponseFunction
217
218   const Int_t kRawSignalOverflow = 0x3FF ; // decimal 1023
219   Bool_t lowGain = kFALSE ; 
220
221   TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
222
223   for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
224     signalF.SetParameter(0, fgHighCharge) ; 
225     signalF.SetParameter(1, fgHighGain) ; 
226     signalF.SetParameter(2, fAmplitude) ; 
227     signalF.SetParameter(3, fTZero) ; 
228     Double_t time = iTime * GetRawFormatTimeMax() / GetRawFormatTimeBins() ;
229     Double_t signal = signalF.Eval(time) ;     
230     fDataHG[iTime] += signal;
231     if ( static_cast<Int_t>(fDataHG[iTime]+0.5) > kRawSignalOverflow ){  // larger than 10 bits 
232       fDataHG[iTime] = kRawSignalOverflow ;
233       lowGain = kTRUE ; 
234     }
235
236     signalF.SetParameter(0, GetRawFormatLowCharge() ) ;     
237     signalF.SetParameter(1, GetRawFormatLowGain() ) ; 
238     signal = signalF.Eval(time) ;  
239     fDataLG[iTime] += signal;
240     if ( static_cast<Int_t>(fDataLG[iTime]+0.5) > kRawSignalOverflow)  // larger than 10 bits 
241       fDataLG[iTime] = kRawSignalOverflow ;
242
243   }
244   // Digitize floating point amplitudes to integers
245   if (fDigitize) Digitize();
246   return lowGain ; 
247 }
248
249 //__________________________________________________________________
250 void AliPHOSPulseGenerator::GetSamples(Int_t *adcHG, Int_t *adcLG) const
251 {
252   // Return integer sample arrays adcHG and adcLG
253   for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
254     adcHG[iTime] = static_cast<Int_t>(fDataHG[iTime]) ;
255     adcLG[iTime] = static_cast<Int_t>(fDataLG[iTime]) ;
256   }
257 }
258
259 //__________________________________________________________________
260 void AliPHOSPulseGenerator::Print(Option_t*) const
261 {
262   // Prints sampled amplitudes to stdout
263   Int_t i;
264   cout << "High gain: ";
265   for (i=0; i<fkTimeBins; i++)
266     cout << (Int_t)fDataHG[i] << " ";
267   cout << endl;
268
269   cout << "Low  gain: ";
270   for (i=0; i<fkTimeBins; i++)
271     cout << (Int_t)fDataLG[i] << " ";
272   cout << endl;
273 }
274
275 //__________________________________________________________________
276 void AliPHOSPulseGenerator::Draw(Option_t*)
277 {
278   // Draw graphs with high and low gain samples
279
280   Double_t *time = new Double_t[fkTimeBins];
281   for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
282     time[iTime] = iTime * GetRawFormatTimeMax() / GetRawFormatTimeBins() ;
283   }
284   Int_t nPoints = fkTimeBins;
285   TGraph *graphHG = new TGraph(nPoints,time,fDataHG);
286   TGraph *graphLG = new TGraph(nPoints,time,fDataLG);
287   graphHG->SetMarkerStyle(20);
288   graphLG->SetMarkerStyle(20);
289   graphHG->SetMarkerSize(0.4);
290   graphLG->SetMarkerSize(0.4);
291   graphHG->SetTitle("High gain samples");
292   graphLG->SetTitle("Low gain samples");
293
294   TCanvas *c1 = new TCanvas("c1","Raw ALTRO samples",10,10,700,500);
295   c1->SetFillColor(0);
296   c1->Divide(2,1);
297   c1->cd(1);
298   gPad->SetLeftMargin(0.15);
299   graphHG->Draw("AP");
300   graphHG->GetHistogram()->SetTitleOffset(1.6,"Y");
301   graphHG->GetHistogram()->SetXTitle("time, #musec");
302   graphHG->GetHistogram()->SetYTitle("Amplitude, ADC counts");
303   c1->cd(2);
304   gPad->SetLeftMargin(0.15);
305   graphLG->Draw("AP");
306   graphLG->GetHistogram()->SetTitleOffset(1.6,"Y");
307   graphLG->GetHistogram()->SetXTitle("time, #musec");
308   graphLG->GetHistogram()->SetYTitle("Amplitude, ADC counts");
309   c1->Update();
310 }