1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 /* History of cvs commits:
20 * Revision 1.49 2007/01/22 17:29:12 pavlinov
21 * EMCAL geometry can be created independently form anything now
23 * Revision 1.48 2006/12/19 02:34:13 pavlinov
24 * clean up the EMCAL name scheme : super module -> module -> tower (or cell)
26 * Revision 1.47 2006/12/05 17:12:03 gustavo
27 * Updated AliEMCAL::Digits2Raw, reads first provisional RCU mapping files to make Raw data with new AliCaloAltroMapping and AliCaloRawStream
31 //_________________________________________________________________________
32 // Base Class for EMCAL description:
33 // This class contains material definitions
34 // for the EMCAL - It does not place the detector in Alice
35 //*-- Author: Yves Schutz (SUBATECH)
37 //*-- Additional Contributions: Sahal Yacoob (LBNL/UCT)
39 //////////////////////////////////////////////////////////////////////////////
41 // --- ROOT system ---
45 #include <TVirtualMC.h>
51 // --- Standard library ---
53 // --- AliRoot header files ---
57 #include "AliEMCALLoader.h"
58 #include "AliEMCALSDigitizer.h"
59 #include "AliEMCALDigitizer.h"
60 #include "AliEMCALDigit.h"
61 //#include "AliAltroMapping.h"
62 #include "AliCaloAltroMapping.h"
63 #include "AliAltroBuffer.h"
64 #include "AliRawReader.h"
65 #include "AliCaloRawStream.h"
69 Double_t AliEMCAL::fgCapa = 1.; // 1pF
70 Int_t AliEMCAL::fgOrder = 2 ;
71 Double_t AliEMCAL::fgTimeMax = 2.56E-5 ; // each sample is over 100 ns fTimeMax/fTimeBins
72 Double_t AliEMCAL::fgTimePeak = 4.1E-6 ; // 4 micro seconds
73 Double_t AliEMCAL::fgTimeTrigger = 100E-9 ; // 100ns, just for a reference
74 // some digitization constants
75 Int_t AliEMCAL::fgThreshold = 1;
76 Int_t AliEMCAL::fgDDLPerSuperModule = 2; // 2 ddls per SuperModule
78 //____________________________________________________________________________
86 fHighLowGainFactor(0.),
95 //____________________________________________________________________________
96 AliEMCAL::AliEMCAL(const char* name, const char* title)
97 : AliDetector(name,title),
103 fHighLowGainFactor(0.),
106 // ctor : title is used to identify the layout
111 //____________________________________________________________________________
112 AliEMCAL::~AliEMCAL()
117 //____________________________________________________________________________
118 void AliEMCAL::InitConstants()
120 //initialize EMCAL values
122 fBirkC1 = 0.013/1.032;
123 fBirkC2 = 9.6e-6/(1.032 * 1.032);
125 fHighCharge = 8.2 ; // adjusted for a high gain range of 5.12 GeV (10 bits)
127 fHighLowGainFactor = 16. ; // adjusted for a low gain range of 82 GeV (10 bits)
128 fLowGainOffset = 1 ; // offset added to the module id to distinguish high and low gain data
131 //____________________________________________________________________________
132 void AliEMCAL::Init()
135 Int_t * idtmed = fIdtmed->GetArray() - 1599 ;
136 // --- Set decent energy thresholds for gamma and electron tracking
138 // Tracking threshold for photons and electrons in Lead
139 Float_t cutgam=10.e-5; // 100 kev;
140 Float_t cutele=10.e-5; // 100 kev;
141 TString ntmp(GetTitle());
143 if(ntmp.Contains("10KEV")) {
144 cutele = cutgam = 1.e-5;
145 } else if(ntmp.Contains("50KEV")) {
146 cutele = cutgam = 5.e-5;
147 } else if(ntmp.Contains("100KEV")) {
148 cutele = cutgam = 1.e-4;
149 } else if(ntmp.Contains("200KEV")) {
150 cutele = cutgam = 2.e-4;
151 } else if(ntmp.Contains("500KEV")) {
152 cutele = cutgam = 5.e-4;
155 gMC->Gstpar(idtmed[1600],"CUTGAM", cutgam);
156 gMC->Gstpar(idtmed[1600],"CUTELE", cutele); // 1MEV -> 0.1MEV; 15-aug-05
157 gMC->Gstpar(idtmed[1600],"BCUTE", cutgam); // BCUTE and BCUTM start from GUTGUM
158 gMC->Gstpar(idtmed[1600],"BCUTM", cutgam); // BCUTE and BCUTM start from GUTGUM
159 // --- Generate explicitly delta rays in Lead ---
160 gMC->Gstpar(idtmed[1600], "LOSS", 3) ;
161 gMC->Gstpar(idtmed[1600], "DRAY", 1) ;
162 gMC->Gstpar(idtmed[1600], "DCUTE", cutele) ;
163 gMC->Gstpar(idtmed[1600], "DCUTM", cutele) ;
165 // --- in aluminium parts ---
166 gMC->Gstpar(idtmed[1602],"CUTGAM", cutgam) ;
167 gMC->Gstpar(idtmed[1602],"CUTELE", cutele) ;
168 gMC->Gstpar(idtmed[1602],"BCUTE", cutgam); // BCUTE and BCUTM start from GUTGUM
169 gMC->Gstpar(idtmed[1602],"BCUTM", cutgam); // BCUTE and BCUTM start from GUTGUM
170 gMC->Gstpar(idtmed[1602], "LOSS",3.) ;
171 gMC->Gstpar(idtmed[1602], "DRAY",1.) ;
172 gMC->Gstpar(idtmed[1602], "DCUTE", cutele) ;
173 gMC->Gstpar(idtmed[1602], "DCUTM", cutele) ;
175 // --- and finally thresholds for photons and electrons in the scintillator ---
176 gMC->Gstpar(idtmed[1601],"CUTGAM", cutgam) ;
177 gMC->Gstpar(idtmed[1601],"CUTELE", cutele) ;// 1MEV -> 0.1MEV; 15-aug-05
178 gMC->Gstpar(idtmed[1601],"BCUTE", cutgam); // BCUTE and BCUTM start from GUTGUM
179 gMC->Gstpar(idtmed[1601],"BCUTM", cutgam); // BCUTE and BCUTM start from GUTGUM
180 gMC->Gstpar(idtmed[1601], "LOSS",3) ; // generate delta rays
181 gMC->Gstpar(idtmed[1601], "DRAY",1) ;
182 gMC->Gstpar(idtmed[1601], "DCUTE", cutele) ;
183 gMC->Gstpar(idtmed[1601], "DCUTM", cutele) ;
186 gMC->Gstpar(idtmed[1603],"CUTGAM", cutgam);
187 gMC->Gstpar(idtmed[1603],"CUTELE", cutele);
188 gMC->Gstpar(idtmed[1603],"BCUTE", cutgam); // BCUTE and BCUTM start from GUTGUM
189 gMC->Gstpar(idtmed[1603],"BCUTM", cutgam); // BCUTE and BCUTM start from GUTGUM
190 // --- Generate explicitly delta rays
191 gMC->Gstpar(idtmed[1603], "LOSS",3);
192 gMC->Gstpar(idtmed[1603], "DRAY",1);
193 gMC->Gstpar(idtmed[1603], "DCUTE", cutele) ;
194 gMC->Gstpar(idtmed[1603], "DCUTM", cutele) ;
196 AliEMCALGeometry* geom = GetGeometry();
197 if(geom->GetILOSS()>=0) {
198 for(int i=1600; i<=1603; i++) gMC->Gstpar(idtmed[i], "LOSS", geom->GetILOSS()) ;
200 if(geom->GetIHADR()>=0) {
201 for(int i=1600; i<=1603; i++) gMC->Gstpar(idtmed[i], "HADR", geom->GetIHADR()) ;
205 //____________________________________________________________________________
206 AliDigitizer* AliEMCAL::CreateDigitizer(AliRunDigitizer* manager) const
208 //create and return the digitizer
209 return new AliEMCALDigitizer(manager);
212 //____________________________________________________________________________
213 void AliEMCAL::CreateMaterials()
215 // Definitions of materials to build EMCAL and associated tracking media.
216 // media number in idtmed are 1599 to 1698.
218 Float_t aAir[4]={12.0107,14.0067,15.9994,39.948};
219 Float_t zAir[4]={6.,7.,8.,18.};
220 Float_t wAir[4]={0.000124,0.755267,0.231781,0.012827};
221 Float_t dAir = 1.20479E-3;
222 AliMixture(0, "Air$", aAir, zAir, dAir, 4, wAir) ;
225 AliMaterial(1, "Pb$", 207.2, 82, 11.35, 0.56, 0., 0, 0) ;
228 // --- The polysterene scintillator (CH) ---
229 Float_t aP[2] = {12.011, 1.00794} ;
230 Float_t zP[2] = {6.0, 1.0} ;
231 Float_t wP[2] = {1.0, 1.0} ;
234 AliMixture(2, "Polystyrene$", aP, zP, dP, -2, wP) ;
237 AliMaterial(3, "Al$", 26.98, 13., 2.7, 8.9, 999., 0, 0) ;
238 // --- Absorption length is ignored ^
240 // 25-aug-04 by PAI - see PMD/AliPMDv0.cxx for STEEL definition
241 Float_t asteel[4] = { 55.847,51.9961,58.6934,28.0855 };
242 Float_t zsteel[4] = { 26.,24.,28.,14. };
243 Float_t wsteel[4] = { .715,.18,.1,.005 };
244 AliMixture(4, "STAINLESS STEEL$", asteel, zsteel, 7.88, 4, wsteel);
246 // DEFINITION OF THE TRACKING MEDIA
248 // for EMCAL: idtmed[1599->1698] equivalent to fIdtmed[0->100]
249 Int_t isxfld = gAlice->Field()->Integ() ;
250 Float_t sxmgmx = gAlice->Field()->Max() ;
252 // Air -> idtmed[1599]
253 AliMedium(0, "Air$", 0, 0,
254 isxfld, sxmgmx, 10.0, 1.0, 0.1, 0.1, 10.0, 0, 0) ;
256 // The Lead -> idtmed[1600]
258 AliMedium(1, "Lead$", 1, 0,
259 isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
261 // The scintillator of the CPV made of Polystyrene scintillator -> idtmed[1601]
262 float deemax = 0.1; // maximum fractional energy loss in one step (0 < DEEMAX â
\89¤ 1);i
263 AliMedium(2, "Scintillator$", 2, 1,
264 isxfld, sxmgmx, 10.0, 0.001, deemax, 0.001, 0.001, 0, 0) ;
266 // Various Aluminium parts made of Al -> idtmed[1602]
267 AliMedium(3, "Al$", 3, 0,
268 isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0) ;
270 // 25-aug-04 by PAI : see PMD/AliPMDv0.cxx for STEEL definition -> idtmed[1603]
271 AliMedium(4, "S steel$", 4, 0,
272 isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0) ;
275 //set constants for Birk's Law implentation
278 fBirkC2 = 9.6e-6/(dP * dP);
281 //____________________________________________________________________________
282 void AliEMCAL::Digits2Raw()
284 // convert digits of the current event to raw data
286 AliEMCALLoader * loader = dynamic_cast<AliEMCALLoader*>(fLoader) ;
289 loader->LoadDigits("EMCAL");
291 TClonesArray* digits = loader->Digits() ;
294 Error("Digits2Raw", "no digits found !");
299 loader->LoadDigitizer();
300 AliEMCALDigitizer * digitizer = dynamic_cast<AliEMCALDigitizer *>(loader->Digitizer()) ;
303 AliEMCALGeometry* geom = GetGeometry();
305 AliError(Form("No geometry found !"));
309 AliAltroBuffer* buffer = NULL;
311 Int_t adcValuesLow[fgkTimeBins];
312 Int_t adcValuesHigh[fgkTimeBins];
314 //Load Mapping RCU files once
315 TString path = gSystem->Getenv("ALICE_ROOT");
316 path += "/EMCAL/mapping/RCU";
317 TString path0 = path+"0.data";//This file will change in future
318 TString path1 = path+"1.data";//This file will change in future
319 AliAltroMapping * mapping[2] ; // For the moment only 2
320 mapping[0] = new AliCaloAltroMapping(path0.Data());
321 mapping[1] = new AliCaloAltroMapping(path1.Data());
323 // loop over digits (assume ordered digits)
324 for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
325 AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
326 if (digit->GetAmp() < fgThreshold)
336 geom->GetCellIndex(digit->GetId(), nSM, nModule, nIphi, nIeta);
337 geom->GetCellPhiEtaIndexInSModule(nSM, nModule, nIphi, nIeta,iphi, ieta) ;
339 //Check which is the RCU of the cell.
342 if (0<=iphi&&iphi<8) iRCU=0; // first cable row
343 else if (8<=iphi&&iphi<16 && 0<=ieta&&ieta<24) iRCU=0; // first half;
346 else if(8<=iphi&&iphi<16 && 24<=ieta&&ieta<48) iRCU=1; // second half;
348 else if(16<=iphi&&iphi<24) iRCU=1; // third cable row
351 Int_t iDDL = fgDDLPerSuperModule* nSM + iRCU;
354 if (iDDL != prevDDL) {
355 // write real header and close previous file
359 buffer->WriteDataHeader(kFALSE, kFALSE);
363 // open new file and write dummy header
364 TString fileName = AliDAQ::DdlFileName("EMCAL",iDDL);
365 buffer = new AliAltroBuffer(fileName.Data(),mapping[iRCU]);
366 buffer->WriteDataHeader(kTRUE, kFALSE); //Dummy;
370 // out of time range signal (?)
371 if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
372 AliInfo("Signal is out of time range.\n");
373 buffer->FillBuffer((Int_t)digit->GetAmp());
374 buffer->FillBuffer(GetRawFormatTimeBins() ); // time bin
375 buffer->FillBuffer(3); // bunch length
376 buffer->WriteTrailer(3, ieta, iphi, nSM); // trailer
377 // calculate the time response function
380 Double_t energy = digit->GetAmp() * digitizer->GetECAchannel() + digitizer->GetECApedestal() ;
382 Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), energy, adcValuesHigh, adcValuesLow) ;
385 buffer->WriteChannel(ieta, iphi, 0, GetRawFormatTimeBins(), adcValuesLow, fgThreshold);
387 buffer->WriteChannel(ieta,iphi, 1, GetRawFormatTimeBins(), adcValuesHigh, fgThreshold);
392 // write real header and close last file
395 buffer->WriteDataHeader(kFALSE, kFALSE);
398 mapping[0]->Delete();
399 mapping[1]->Delete();
400 loader->UnloadDigits();
403 //____________________________________________________________________________
404 void AliEMCAL::Raw2Digits(AliRawReader* reader)
406 // convert raw data of the current event to digits
407 AliEMCALGeometry * geom = GetGeometry();
408 AliEMCALLoader * loader = dynamic_cast<AliEMCALLoader*>(fLoader) ;
411 loader->CleanDigits(); // start from scratch
412 loader->LoadDigits("EMCAL");
413 TClonesArray* digits = loader->Digits() ;
414 digits->Clear(); // yes, this is perhaps somewhat paranoid.. [clearing an extra time]
417 Error("Raw2Digits", "no digits found !");
421 Error("Raw2Digits", "no raw reader found !");
425 // and get the digitizer too
426 loader->LoadDigitizer();
427 AliEMCALDigitizer * digitizer = dynamic_cast<AliEMCALDigitizer *>(loader->Digitizer()) ;
429 // Use AliAltroRawStream to read the ALTRO format. No need to
430 // reinvent the wheel :-)
431 AliCaloRawStream in(reader,"EMCAL");
432 // Select EMCAL DDL's;
433 reader->Select("EMCAL");
435 // reading is from previously existing AliEMCALGetter.cxx
437 Bool_t first = kTRUE ;
439 TF1 * signalF = new TF1("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
440 signalF->SetParNames("Charge", "Gain", "Amplitude", "TimeZero");
442 Bool_t lowGainFlag = kFALSE ;
447 Double_t energy = 0. ;
449 TGraph * gLowGain = new TGraph(GetRawFormatTimeBins()) ;
450 TGraph * gHighGain= new TGraph(GetRawFormatTimeBins()) ;
452 while ( in.Next() ) { // EMCAL entries loop
453 if ( in.IsNewRow() ) {//phi
454 if ( in.IsNewColumn() ) {//eta
455 id = geom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
457 FitRaw(lowGainFlag, gLowGain, gHighGain, signalF, energy, time) ;
459 if (time == 0. && energy == 0.) {
463 amp = static_cast<Int_t>( (energy - digitizer->GetECApedestal()) / digitizer->GetECAchannel() + 0.5 ) ;
467 new((*digits)[idigit]) AliEMCALDigit( -1, -1, id, amp, time) ;
471 for (index = 0; index < GetRawFormatTimeBins(); index++) {
472 gLowGain->SetPoint(index, index * GetRawFormatTimeMax() / GetRawFormatTimeBins(), 0) ;
473 gHighGain->SetPoint(index, index * GetRawFormatTimeMax() / GetRawFormatTimeBins(), 0) ;
477 id = geom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
478 if (in.GetModule() == GetRawFormatLowGainOffset() ) {
479 lowGainFlag = kTRUE ;
482 lowGainFlag = kFALSE ;
487 gLowGain->SetPoint(in.GetTime(),
488 in.GetTime()* GetRawFormatTimeMax() / GetRawFormatTimeBins(),
492 gHighGain->SetPoint(in.GetTime(),
493 in.GetTime() * GetRawFormatTimeMax() / GetRawFormatTimeBins(),
496 } // EMCAL entries loop
506 //____________________________________________________________________________
507 void AliEMCAL::FitRaw(Bool_t lowGainFlag, TGraph * gLowGain, TGraph * gHighGain, TF1* signalF, Double_t & energy, Double_t & time)
509 // Fits the raw signal time distribution; from AliEMCALGetter
511 const Int_t kNoiseThreshold = 0 ;
512 Double_t timezero1 = 0., timezero2 = 0., timemax = 0. ;
513 Double_t signal = 0., signalmax = 0. ;
517 timezero1 = timezero2 = signalmax = timemax = 0. ;
518 signalF->FixParameter(0, GetRawFormatLowCharge()) ;
519 signalF->FixParameter(1, GetRawFormatLowGain()) ;
521 for (index = 0; index < GetRawFormatTimeBins(); index++) {
522 gLowGain->GetPoint(index, time, signal) ;
523 if (signal > kNoiseThreshold && timezero1 == 0.)
525 if (signal <= kNoiseThreshold && timezero1 > 0. && timezero2 == 0.)
527 if (signal > signalmax) {
532 signalmax /= RawResponseFunctionMax(GetRawFormatLowCharge(),
533 GetRawFormatLowGain()) ;
534 if ( timezero1 + GetRawFormatTimePeak() < GetRawFormatTimeMax() * 0.4 ) { // else its noise
535 signalF->SetParameter(2, signalmax) ;
536 signalF->SetParameter(3, timezero1) ;
537 gLowGain->Fit(signalF, "QRON", "", 0., timezero2); //, "QRON") ;
538 energy = signalF->GetParameter(2) ;
539 time = signalF->GetMaximumX() - GetRawFormatTimePeak() - GetRawFormatTimeTrigger() ;
542 timezero1 = timezero2 = signalmax = timemax = 0. ;
543 signalF->FixParameter(0, GetRawFormatHighCharge()) ;
544 signalF->FixParameter(1, GetRawFormatHighGain()) ;
546 for (index = 0; index < GetRawFormatTimeBins(); index++) {
547 gHighGain->GetPoint(index, time, signal) ;
548 if (signal > kNoiseThreshold && timezero1 == 0.)
550 if (signal <= kNoiseThreshold && timezero1 > 0. && timezero2 == 0.)
552 if (signal > signalmax) {
557 signalmax /= RawResponseFunctionMax(GetRawFormatHighCharge(),
558 GetRawFormatHighGain()) ;;
559 if ( timezero1 + GetRawFormatTimePeak() < GetRawFormatTimeMax() * 0.4 ) { // else its noise
560 signalF->SetParameter(2, signalmax) ;
561 signalF->SetParameter(3, timezero1) ;
562 gHighGain->Fit(signalF, "QRON", "", 0., timezero2) ;
563 energy = signalF->GetParameter(2) ;
564 time = signalF->GetMaximumX() - GetRawFormatTimePeak() - GetRawFormatTimeTrigger() ;
571 //____________________________________________________________________________
572 void AliEMCAL::Hits2SDigits()
574 // create summable digits
577 AliEMCALSDigitizer emcalDigitizer(fLoader->GetRunLoader()->GetFileName().Data()) ;
578 emcalDigitizer.SetEventRange(0, -1) ; // do all the events
579 emcalDigitizer.ExecuteTask() ;
582 //____________________________________________________________________________
584 AliLoader* AliEMCAL::MakeLoader(const char* topfoldername)
586 //different behaviour than standard (singleton getter)
587 // --> to be discussed and made eventually coherent
588 fLoader = new AliEMCALLoader(GetName(),topfoldername);
592 //__________________________________________________________________
593 Double_t AliEMCAL::RawResponseFunction(Double_t *x, Double_t *par)
595 // Shape of the electronics raw reponse:
596 // It is a semi-gaussian, 2nd order Gamma function of the general form
597 // v(t) = n**n * Q * A**n / C *(t/tp)**n * exp(-n * t/tp) with
598 // tp : peaking time par[0]
599 // n : order of the function
600 // C : integrating capacitor in the preamplifier
601 // A : open loop gain of the preamplifier
602 // Q : the total APD charge to be measured Q = C * energy
605 Double_t xx = x[0] - ( fgTimeTrigger + par[3] ) ;
607 if (xx < 0 || xx > fgTimeMax)
610 Double_t fac = par[0] * TMath::Power(fgOrder, fgOrder) * TMath::Power(par[1], fgOrder) / fgCapa ;
611 signal = fac * par[2] * TMath::Power(xx / fgTimePeak, fgOrder) * TMath::Exp(-fgOrder * (xx / fgTimePeak)) ;
616 //__________________________________________________________________
617 Double_t AliEMCAL::RawResponseFunctionMax(Double_t charge, Double_t gain)
619 //compute the maximum of the raw response function and return
620 return ( charge * TMath::Power(fgOrder, fgOrder) * TMath::Power(gain, fgOrder)
621 / ( fgCapa * TMath::Exp(fgOrder) ) );
624 //__________________________________________________________________
625 Bool_t AliEMCAL::RawSampledResponse(
626 const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL) const
628 // for a start time dtime and an amplitude damp given by digit,
629 // calculates the raw sampled response AliEMCAL::RawResponseFunction
631 const Int_t kRawSignalOverflow = 0x3FF ;
632 Bool_t lowGain = kFALSE ;
634 TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
636 for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
637 signalF.SetParameter(0, GetRawFormatHighCharge() ) ;
638 signalF.SetParameter(1, GetRawFormatHighGain() ) ;
639 signalF.SetParameter(2, damp) ;
640 signalF.SetParameter(3, dtime) ;
641 Double_t time = iTime * GetRawFormatTimeMax() / GetRawFormatTimeBins() ;
642 Double_t signal = signalF.Eval(time) ;
643 if ( static_cast<Int_t>(signal+0.5) > kRawSignalOverflow ){ // larger than 10 bits
644 signal = kRawSignalOverflow ;
647 adcH[iTime] = static_cast<Int_t>(signal + 0.5) ;
649 signalF.SetParameter(0, GetRawFormatLowCharge() ) ;
650 signalF.SetParameter(1, GetRawFormatLowGain() ) ;
651 signal = signalF.Eval(time) ;
652 if ( static_cast<Int_t>(signal+0.5) > kRawSignalOverflow) // larger than 10 bits
653 signal = kRawSignalOverflow ;
654 adcL[iTime] = static_cast<Int_t>(0.5 + signal ) ;