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