]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EMCAL/AliEMCALRawUtils.cxx
Linkning with libMinuit on Windows/Cygwin
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALRawUtils.cxx
CommitLineData
ee299369 1/**************************************************************************
2 * Copyright(c) 1998-1999, 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/* History of cvs commits:
18 *
c47157cd 19 * $Log$
48a56166 20 * Revision 1.10 2007/12/06 13:58:11 hristov
21 * Additional pritection. Do not delete the mapping, it is owned by another class
22 *
e36e3bcf 23 * Revision 1.9 2007/12/06 02:19:51 jklay
24 * incorporated fitting procedure from testbeam analysis into AliRoot
25 *
8cb998bd 26 * Revision 1.8 2007/12/05 02:30:51 jklay
27 * modification to read Altro mappings into AliEMCALRecParam and pass to AliEMCALRawUtils from AliEMCALReconstructor; add option to AliEMCALRawUtils to set old RCU format (for testbeam) or not
28 *
feedcab9 29 * Revision 1.7 2007/11/14 15:51:46 gustavo
30 * Take out few unnecessary prints
31 *
6d0b6861 32 * Revision 1.6 2007/11/01 01:23:51 mvl
33 * Removed call to SetOldRCUFormat, which is only needed for testbeam data
34 *
f31c5945 35 * Revision 1.5 2007/11/01 01:20:33 mvl
36 * Further improvement of peak finding; more robust fit
37 *
e9dbb64a 38 * Revision 1.4 2007/10/31 17:15:24 mvl
39 * Fixed bug in raw data unpacking; Added pedestal to signal fit; Added logic to deal with high/low gain
40 *
82cbdfca 41 * Revision 1.3 2007/09/27 08:36:46 mvl
42 * More robust setting of fit range in FitRawSignal (P. Hristov)
43 *
5a056daa 44 * Revision 1.2 2007/09/03 20:55:35 jklay
45 * EMCAL e-by-e reconstruction methods from Cvetan
46 *
c47157cd 47 * Revision 1.1 2007/03/17 19:56:38 mvl
48 * Moved signal shape routines from AliEMCAL to separate class AliEMCALRawUtils to streamline raw data reconstruction code.
49 * */
ee299369 50
51//*-- Author: Marco van Leeuwen (LBL)
52#include "AliEMCALRawUtils.h"
53
54#include "TF1.h"
55#include "TGraph.h"
56#include "TSystem.h"
57
58#include "AliLog.h"
72c58de0 59#include "AliRun.h"
ee299369 60#include "AliRunLoader.h"
61#include "AliCaloAltroMapping.h"
62#include "AliAltroBuffer.h"
63#include "AliRawReader.h"
64#include "AliCaloRawStream.h"
65#include "AliDAQ.h"
66
feedcab9 67#include "AliEMCALRecParam.h"
ee299369 68#include "AliEMCALLoader.h"
69#include "AliEMCALGeometry.h"
70#include "AliEMCALDigitizer.h"
71#include "AliEMCALDigit.h"
20b636fc 72#include "AliEMCAL.h"
ee299369 73
74ClassImp(AliEMCALRawUtils)
75
76// Signal shape parameters
82cbdfca 77Double_t AliEMCALRawUtils::fgTimeBinWidth = 100E-9 ; // each sample is 100 ns
ee299369 78Double_t AliEMCALRawUtils::fgTimeTrigger = 1.5E-6 ; // 15 time bins ~ 1.5 musec
79
80// some digitization constants
81Int_t AliEMCALRawUtils::fgThreshold = 1;
82Int_t AliEMCALRawUtils::fgDDLPerSuperModule = 2; // 2 ddls per SuperModule
fe93d365 83Int_t AliEMCALRawUtils::fgPedestalValue = 32; // pedestal value for digits2raw
84Double_t AliEMCALRawUtils::fgFEENoise = 3.; // 3 ADC channels of noise (sampled)
ee299369 85
8cb998bd 86AliEMCALRawUtils::AliEMCALRawUtils()
b4133f05 87 : fHighLowGainFactor(0.), fOrder(0), fTau(0.), fNoiseThreshold(0),
88 fNPedSamples(0), fGeom(0), fOption("")
8cb998bd 89{
b4133f05 90
91 //These are default parameters.
92 //Can be re-set from without with setter functions
ee299369 93 fHighLowGainFactor = 16. ; // adjusted for a low gain range of 82 GeV (10 bits)
b4133f05 94 fOrder = 2; // order of gamma fn
95 fTau = 2.35; // in units of timebin, from CERN 2007 testbeam
96 fNoiseThreshold = 3;
97 fNPedSamples = 5;
65bdc82f 98
99 //Get Mapping RCU files from the AliEMCALRecParam
100 const TObjArray* maps = AliEMCALRecParam::GetMappings();
101 if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
102
103 for(Int_t i = 0; i < 2; i++) {
104 fMapping[i] = (AliAltroMapping*)maps->At(i);
105 }
106
72c58de0 107 //To make sure we match with the geometry in a simulation file,
108 //let's try to get it first. If not, take the default geometry
109 AliRunLoader *rl = AliRunLoader::GetRunLoader();
110 if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL")) {
111 fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
112 } else {
113 AliInfo(Form("Using default geometry in raw reco"));
937d0661 114 fGeom = AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
65bdc82f 115 }
116
72c58de0 117 if(!fGeom) AliFatal(Form("Could not get geometry!"));
118
65bdc82f 119}
120
121//____________________________________________________________________________
122AliEMCALRawUtils::AliEMCALRawUtils(const AliEMCALRawUtils& rawU)
123 : TObject(),
124 fHighLowGainFactor(rawU.fHighLowGainFactor),
b4133f05 125 fOrder(rawU.fOrder),
126 fTau(rawU.fTau),
127 fNoiseThreshold(rawU.fNoiseThreshold),
128 fNPedSamples(rawU.fNPedSamples),
65bdc82f 129 fGeom(rawU.fGeom),
130 fOption(rawU.fOption)
131{
132 //copy ctor
133 fMapping[0] = rawU.fMapping[0];
134 fMapping[1] = rawU.fMapping[1];
135}
136
137//____________________________________________________________________________
138AliEMCALRawUtils& AliEMCALRawUtils::operator =(const AliEMCALRawUtils &rawU)
139{
140 //assignment operator
141
142 if(this != &rawU) {
143 fHighLowGainFactor = rawU.fHighLowGainFactor;
b4133f05 144 fOrder = rawU.fOrder;
145 fTau = rawU.fTau;
146 fNoiseThreshold = rawU.fNoiseThreshold;
147 fNPedSamples = rawU.fNPedSamples;
65bdc82f 148 fGeom = rawU.fGeom;
149 fOption = rawU.fOption;
150 fMapping[0] = rawU.fMapping[0];
151 fMapping[1] = rawU.fMapping[1];
152 }
153
154 return *this;
155
ee299369 156}
65bdc82f 157
ee299369 158//____________________________________________________________________________
159AliEMCALRawUtils::~AliEMCALRawUtils() {
65bdc82f 160
ee299369 161}
65bdc82f 162
ee299369 163//____________________________________________________________________________
65bdc82f 164void AliEMCALRawUtils::Digits2Raw()
ee299369 165{
166 // convert digits of the current event to raw data
167
168 AliRunLoader *rl = AliRunLoader::GetRunLoader();
169 AliEMCALLoader *loader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
170
171 // get the digits
172 loader->LoadDigits("EMCAL");
173 loader->GetEvent();
174 TClonesArray* digits = loader->Digits() ;
175
176 if (!digits) {
177 Warning("Digits2Raw", "no digits found !");
178 return;
179 }
65bdc82f 180
ee299369 181 static const Int_t nDDL = 12*2; // 12 SM hardcoded for now. Buffers allocated dynamically, when needed, so just need an upper limit here
182 AliAltroBuffer* buffers[nDDL];
183 for (Int_t i=0; i < nDDL; i++)
184 buffers[i] = 0;
185
186 Int_t adcValuesLow[fgkTimeBins];
187 Int_t adcValuesHigh[fgkTimeBins];
188
ee299369 189 // loop over digits (assume ordered digits)
190 for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
191 AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
192 if (digit->GetAmp() < fgThreshold)
193 continue;
194
195 //get cell indices
196 Int_t nSM = 0;
197 Int_t nIphi = 0;
198 Int_t nIeta = 0;
199 Int_t iphi = 0;
200 Int_t ieta = 0;
201 Int_t nModule = 0;
65bdc82f 202 fGeom->GetCellIndex(digit->GetId(), nSM, nModule, nIphi, nIeta);
203 fGeom->GetCellPhiEtaIndexInSModule(nSM, nModule, nIphi, nIeta,iphi, ieta) ;
ee299369 204
205 //Check which is the RCU of the cell.
206 Int_t iRCU = -111;
207 //RCU0
208 if (0<=iphi&&iphi<8) iRCU=0; // first cable row
209 else if (8<=iphi&&iphi<16 && 0<=ieta&&ieta<24) iRCU=0; // first half;
210 //second cable row
211 //RCU1
212 else if(8<=iphi&&iphi<16 && 24<=ieta&&ieta<48) iRCU=1; // second half;
213 //second cable row
214 else if(16<=iphi&&iphi<24) iRCU=1; // third cable row
e36e3bcf 215 if (iRCU<0)
216 Fatal("Digits2Raw()","Non-existent RCU number: %d", iRCU);
ee299369 217
218 //Which DDL?
219 Int_t iDDL = fgDDLPerSuperModule* nSM + iRCU;
220 if (iDDL >= nDDL)
221 Fatal("Digits2Raw()","Non-existent DDL board number: %d", iDDL);
222
223 if (buffers[iDDL] == 0) {
224 // open new file and write dummy header
225 TString fileName = AliDAQ::DdlFileName("EMCAL",iDDL);
65bdc82f 226 buffers[iDDL] = new AliAltroBuffer(fileName.Data(),fMapping[iRCU]);
ee299369 227 buffers[iDDL]->WriteDataHeader(kTRUE, kFALSE); //Dummy;
228 }
229
230 // out of time range signal (?)
231 if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
232 AliInfo("Signal is out of time range.\n");
233 buffers[iDDL]->FillBuffer((Int_t)digit->GetAmp());
234 buffers[iDDL]->FillBuffer(GetRawFormatTimeBins() ); // time bin
235 buffers[iDDL]->FillBuffer(3); // bunch length
236 buffers[iDDL]->WriteTrailer(3, ieta, iphi, nSM); // trailer
237 // calculate the time response function
238 } else {
239 Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), digit->GetAmp(), adcValuesHigh, adcValuesLow) ;
240 if (lowgain)
241 buffers[iDDL]->WriteChannel(ieta, iphi, 0, GetRawFormatTimeBins(), adcValuesLow, fgThreshold);
242 else
243 buffers[iDDL]->WriteChannel(ieta,iphi, 1, GetRawFormatTimeBins(), adcValuesHigh, fgThreshold);
244 }
245 }
246
247 // write headers and close files
248 for (Int_t i=0; i < nDDL; i++) {
249 if (buffers[i]) {
250 buffers[i]->Flush();
251 buffers[i]->WriteDataHeader(kFALSE, kFALSE);
252 delete buffers[i];
253 }
254 }
65bdc82f 255
ee299369 256 loader->UnloadDigits();
257}
258
259//____________________________________________________________________________
65bdc82f 260void AliEMCALRawUtils::Raw2Digits(AliRawReader* reader,TClonesArray *digitsArr)
ee299369 261{
65bdc82f 262 // convert raw data of the current event to digits
ee299369 263
c47157cd 264 digitsArr->Clear();
ee299369 265
c47157cd 266 if (!digitsArr) {
ee299369 267 Error("Raw2Digits", "no digits found !");
268 return;
269 }
270 if (!reader) {
271 Error("Raw2Digits", "no raw reader found !");
272 return;
273 }
274
65bdc82f 275 AliCaloRawStream in(reader,"EMCAL",fMapping);
ee299369 276 // Select EMCAL DDL's;
277 reader->Select("EMCAL");
feedcab9 278
8cb998bd 279 //Updated fitting routine from 2007 beam test takes into account
280 //possibility of two peaks in data and selects first one for fitting
281 //Also sets some of the starting parameters based on the shape of the
282 //given raw signal being fit
283
284 TF1 * signalF = new TF1("signal", RawResponseFunction, 0, GetRawFormatTimeBins(), 5);
b4133f05 285 signalF->SetParameters(10.,0.,fTau,fOrder,5.); //set all defaults once, just to be safe
8cb998bd 286 signalF->SetParNames("amp","t0","tau","N","ped");
b4133f05 287 signalF->SetParameter(2,fTau); // tau in units of time bin
8cb998bd 288 signalF->SetParLimits(2,2,-1);
b4133f05 289 signalF->SetParameter(3,fOrder); // order
8cb998bd 290 signalF->SetParLimits(3,2,-1);
48a56166 291
ee299369 292 Int_t id = -1;
82cbdfca 293 Float_t time = 0. ;
294 Float_t amp = 0. ;
ee299369 295
8cb998bd 296 //Graph to hold data we will fit (should be converted to an array
297 //later to speed up processing
298 TGraph * gSig = new TGraph(GetRawFormatTimeBins());
ee299369 299
82cbdfca 300 Int_t readOk = 1;
ee299369 301 Int_t lowGain = 0;
302
82cbdfca 303 while (readOk && in.GetModule() < 0)
304 readOk = in.Next(); // Go to first digit
e9dbb64a 305
8cb998bd 306 Int_t col = 0;
307 Int_t row = 0;
308
82cbdfca 309 while (readOk) {
65bdc82f 310 id = fGeom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
ee299369 311 lowGain = in.IsLowGain();
82cbdfca 312 Int_t maxTime = in.GetTime(); // timebins come in reverse order
e9dbb64a 313 if (maxTime < 0 || maxTime >= GetRawFormatTimeBins()) {
314 AliWarning(Form("Invalid time bin %d",maxTime));
315 maxTime = GetRawFormatTimeBins();
316 }
82cbdfca 317 gSig->Set(maxTime+1);
318 // There is some kind of zero-suppression in the raw data,
319 // so set up the TGraph in advance
320 for (Int_t i=0; i < maxTime; i++) {
8cb998bd 321 gSig->SetPoint(i, i , 0);
82cbdfca 322 }
ee299369 323
82cbdfca 324 Int_t iTime = 0;
ee299369 325 do {
82cbdfca 326 if (in.GetTime() >= gSig->GetN()) {
327 AliWarning("Too many time bins");
328 gSig->Set(in.GetTime());
ee299369 329 }
8cb998bd 330 col = in.GetColumn();
331 row = in.GetRow();
332
333 gSig->SetPoint(in.GetTime(), in.GetTime(), in.GetSignal()) ;
fb070798 334
82cbdfca 335 if (in.GetTime() > maxTime)
336 maxTime = in.GetTime();
ee299369 337 iTime++;
82cbdfca 338 } while ((readOk = in.Next()) && !in.IsNewHWAddress());
ee299369 339
340 FitRaw(gSig, signalF, amp, time) ;
ee299369 341
3ced962c 342 if (amp > 0 && amp < 10000) { //check both high and low end of
343 //result, 10000 is somewhat arbitrary
82cbdfca 344 AliDebug(2,Form("id %d lowGain %d amp %g", id, lowGain, amp));
8cb998bd 345 //cout << "col " << col-40 << " row " << row-8 << " lowGain " << lowGain << " amp " << amp << endl;
fb070798 346
82cbdfca 347 AddDigit(digitsArr, id, lowGain, (Int_t)amp, time);
ee299369 348 }
ee299369 349
350 // Reset graph
82cbdfca 351 for (Int_t index = 0; index < gSig->GetN(); index++) {
8cb998bd 352 gSig->SetPoint(index, index, 0) ;
ee299369 353 }
48a56166 354 // Reset starting parameters for fit function
b4133f05 355 signalF->SetParameters(10.,0.,fTau,fOrder,5.); //reset all defaults just to be safe
48a56166 356
82cbdfca 357 }; // EMCAL entries loop
ee299369 358
359 delete signalF ;
360 delete gSig;
361
362 return ;
363}
364
365//____________________________________________________________________________
82cbdfca 366void AliEMCALRawUtils::AddDigit(TClonesArray *digitsArr, Int_t id, Int_t lowGain, Int_t amp, Float_t time) {
367 //
368 // Add a new digit.
369 // This routine checks whether a digit exists already for this tower
370 // and then decides whether to use the high or low gain info
371 //
372 // Called by Raw2Digits
373
374 AliEMCALDigit *digit = 0, *tmpdigit = 0;
375
376 TIter nextdigit(digitsArr);
377 while (digit == 0 && (tmpdigit = (AliEMCALDigit*) nextdigit())) {
378 if (tmpdigit->GetId() == id)
379 digit = tmpdigit;
380 }
381
382 if (!digit) { // no digit existed for this tower; create one
383 if (lowGain)
384 amp = Int_t(fHighLowGainFactor * amp);
385 Int_t idigit = digitsArr->GetEntries();
386 new((*digitsArr)[idigit]) AliEMCALDigit( -1, -1, id, amp, time, idigit) ;
387 }
388 else { // a digit already exists, check range
b4133f05 389 // (use high gain if signal < cut value, otherwise low gain)
82cbdfca 390 if (lowGain) { // new digit is low gain
b4133f05 391 if (digit->GetAmp() > fgkOverflowCut) { // use if stored digit is out of range
82cbdfca 392 digit->SetAmp(Int_t(fHighLowGainFactor * amp));
393 digit->SetTime(time);
394 }
395 }
b4133f05 396 else if (amp < fgkOverflowCut) { // new digit is high gain; use if not out of range
82cbdfca 397 digit->SetAmp(amp);
398 digit->SetTime(time);
399 }
400 }
401}
402
403//____________________________________________________________________________
404void AliEMCALRawUtils::FitRaw(TGraph * gSig, TF1* signalF, Float_t & amp, Float_t & time)
ee299369 405{
406 // Fits the raw signal time distribution; from AliEMCALGetter
407
ee299369 408 amp = time = 0. ;
82cbdfca 409 Double_t ped = 0;
410 Int_t nPed = 0;
ee299369 411
b4133f05 412 for (Int_t index = 0; index < fNPedSamples; index++) {
e9dbb64a 413 Double_t ttime, signal;
82cbdfca 414 gSig->GetPoint(index, ttime, signal) ;
415 if (signal > 0) {
416 ped += signal;
417 nPed++;
ee299369 418 }
419 }
e9dbb64a 420
82cbdfca 421 if (nPed > 0)
422 ped /= nPed;
e9dbb64a 423 else {
8cb998bd 424 AliWarning("Could not determine pedestal");
82cbdfca 425 ped = 10; // put some small value as first guess
e9dbb64a 426 }
82cbdfca 427
e9dbb64a 428 Int_t max_found = 0;
429 Int_t i_max = 0;
430 Float_t max = -1;
8cb998bd 431 Float_t max_fit = gSig->GetN();
e9dbb64a 432 Float_t min_after_sig = 9999;
8cb998bd 433 Int_t tmin_after_sig = gSig->GetN();
e9dbb64a 434 Int_t n_ped_after_sig = 0;
fb070798 435 Int_t plateau_width = 0;
436 Int_t plateau_start = 9999;
e9dbb64a 437
b4133f05 438 for (Int_t i=fNPedSamples; i < gSig->GetN(); i++) {
e9dbb64a 439 Double_t ttime, signal;
440 gSig->GetPoint(i, ttime, signal) ;
441 if (!max_found && signal > max) {
442 i_max = i;
e9dbb64a 443 max = signal;
444 }
b4133f05 445 else if ( max > ped + fNoiseThreshold ) {
e9dbb64a 446 max_found = 1;
447 min_after_sig = signal;
8cb998bd 448 tmin_after_sig = i;
e9dbb64a 449 }
450 if (max_found) {
451 if ( signal < min_after_sig) {
452 min_after_sig = signal;
8cb998bd 453 tmin_after_sig = i;
e9dbb64a 454 }
455 if (i > tmin_after_sig + 5) { // Two close peaks; end fit at minimum
456 max_fit = tmin_after_sig;
457 break;
458 }
b4133f05 459 if ( signal < ped + fNoiseThreshold)
e9dbb64a 460 n_ped_after_sig++;
461 if (n_ped_after_sig >= 5) { // include 5 pedestal bins after peak
8cb998bd 462 max_fit = i;
e9dbb64a 463 break;
464 }
5a056daa 465 }
fb070798 466 //Add check on plateau
467 if (signal >= fgkRawSignalOverflow - fNoiseThreshold) {
468 if(plateau_width == 0) plateau_start = i;
469 plateau_width++;
470 }
471 }
472
473 if(plateau_width > 0) {
474 for(int j = 0; j < plateau_width-2; j++) {
475 //Note, have to remove the same point N times because after each
476 //remove, the positions of all subsequent points have shifted down
477 //We leave the first and last as anchor points
478 gSig->RemovePoint(plateau_start+1);
479 }
82cbdfca 480 }
5a056daa 481
b4133f05 482 if ( max - ped > fNoiseThreshold ) { // else its noise
e9dbb64a 483 AliDebug(2,Form("Fitting max %d ped %d", max, ped));
8cb998bd 484 signalF->SetRange(0,max_fit);
485
486 if(max-ped > 50)
487 signalF->SetParLimits(2,1,3);
488
489 signalF->SetParameter(4, ped) ;
490 signalF->SetParameter(1, i_max);
491 signalF->SetParameter(0, max);
492
493 gSig->Fit(signalF, "QROW"); // Note option 'W': equal errors on all points
494 amp = signalF->GetParameter(0);
495 time = signalF->GetParameter(1)*GetRawFormatTimeBinWidth() - fgTimeTrigger;
ee299369 496 }
497 return;
498}
499//__________________________________________________________________
500Double_t AliEMCALRawUtils::RawResponseFunction(Double_t *x, Double_t *par)
501{
8cb998bd 502 // Matches version used in 2007 beam test
503 //
ee299369 504 // Shape of the electronics raw reponse:
505 // It is a semi-gaussian, 2nd order Gamma function of the general form
506 //
507 // t' = (t - t0 + tau) / tau
508 // F = A * t**N * exp( N * ( 1 - t) ) for t >= 0
509 // F = 0 for t < 0
510 //
511 // parameters:
8cb998bd 512 // A: par[0] // Amplitude = peak value
513 // t0: par[1]
514 // tau: par[2]
515 // N: par[3]
516 // ped: par[4]
ee299369 517 //
518 Double_t signal ;
8cb998bd 519 Double_t tau =par[2];
520 Double_t N =par[3];
521 Double_t ped = par[4];
522 Double_t xx = ( x[0] - par[1] + tau ) / tau ;
ee299369 523
5a056daa 524 if (xx <= 0)
8cb998bd 525 signal = ped ;
ee299369 526 else {
8cb998bd 527 signal = ped + par[0] * TMath::Power(xx , N) * TMath::Exp(N * (1 - xx )) ;
ee299369 528 }
529 return signal ;
530}
531
532//__________________________________________________________________
533Bool_t AliEMCALRawUtils::RawSampledResponse(
534const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL) const
535{
536 // for a start time dtime and an amplitude damp given by digit,
537 // calculates the raw sampled response AliEMCAL::RawResponseFunction
538
ee299369 539 Bool_t lowGain = kFALSE ;
540
48a56166 541 // A: par[0] // Amplitude = peak value
542 // t0: par[1]
543 // tau: par[2]
544 // N: par[3]
545 // ped: par[4]
546
547 TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 5);
548 signalF.SetParameter(0, damp) ;
549 signalF.SetParameter(1, dtime + fgTimeTrigger) ;
b4133f05 550 signalF.SetParameter(2, fTau) ;
551 signalF.SetParameter(3, fOrder);
fe93d365 552 signalF.SetParameter(4, fgPedestalValue);
ee299369 553
554 for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
82cbdfca 555 Double_t time = iTime * GetRawFormatTimeBinWidth() ;
ee299369 556 Double_t signal = signalF.Eval(time) ;
fe93d365 557
558 //According to Terry Awes, 13-Apr-2008
559 //add gaussian noise in quadrature to each sample
560 //Double_t noise = gRandom->Gaus(0.,fgFEENoise);
561 //signal = sqrt(signal*signal + noise*noise);
562
ee299369 563 adcH[iTime] = static_cast<Int_t>(signal + 0.5) ;
b4133f05 564 if ( adcH[iTime] > fgkRawSignalOverflow ){ // larger than 10 bits
565 adcH[iTime] = fgkRawSignalOverflow ;
ee299369 566 lowGain = kTRUE ;
567 }
568
569 signal /= fHighLowGainFactor;
570
571 adcL[iTime] = static_cast<Int_t>(signal + 0.5) ;
b4133f05 572 if ( adcL[iTime] > fgkRawSignalOverflow) // larger than 10 bits
573 adcL[iTime] = fgkRawSignalOverflow ;
ee299369 574 }
575 return lowGain ;
576}