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