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