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.50 2007/02/05 10:43:25 hristov
21 * Changes for correct initialization of Geant4 (Mihaela)
23 * Revision 1.49 2007/01/22 17:29:12 pavlinov
24 * EMCAL geometry can be created independently form anything now
26 * Revision 1.48 2006/12/19 02:34:13 pavlinov
27 * clean up the EMCAL name scheme : super module -> module -> tower (or cell)
29 * Revision 1.47 2006/12/05 17:12:03 gustavo
30 * Updated AliEMCAL::Digits2Raw, reads first provisional RCU mapping files to make Raw data with new AliCaloAltroMapping and AliCaloRawStream
34 //_________________________________________________________________________
35 // Base Class for EMCAL description:
36 // This class contains material definitions
37 // for the EMCAL - It does not place the detector in Alice
38 //*-- Author: Yves Schutz (SUBATECH)
40 //*-- Additional Contributions: Sahal Yacoob (LBNL/UCT)
42 //////////////////////////////////////////////////////////////////////////////
44 // --- ROOT system ---
48 #include <TVirtualMC.h>
54 // --- Standard library ---
56 // --- AliRoot header files ---
60 #include "AliEMCALLoader.h"
61 #include "AliEMCALSDigitizer.h"
62 #include "AliEMCALDigitizer.h"
63 #include "AliEMCALDigit.h"
64 //#include "AliAltroMapping.h"
65 #include "AliCaloAltroMapping.h"
66 #include "AliAltroBuffer.h"
67 #include "AliRawReader.h"
68 #include "AliCaloRawStream.h"
72 Double_t AliEMCAL::fgCapa = 1.; // 1pF
73 Int_t AliEMCAL::fgOrder = 2 ;
74 Double_t AliEMCAL::fgTimeMax = 2.56E-5 ; // each sample is over 100 ns fTimeMax/fTimeBins
75 Double_t AliEMCAL::fgTimePeak = 4.1E-6 ; // 4 micro seconds
76 Double_t AliEMCAL::fgTimeTrigger = 100E-9 ; // 100ns, just for a reference
77 // some digitization constants
78 Int_t AliEMCAL::fgThreshold = 1;
79 Int_t AliEMCAL::fgDDLPerSuperModule = 2; // 2 ddls per SuperModule
81 //____________________________________________________________________________
89 fHighLowGainFactor(0.),
97 //____________________________________________________________________________
98 AliEMCAL::AliEMCAL(const char* name, const char* title)
99 : AliDetector(name,title),
105 fHighLowGainFactor(0.),
108 // ctor : title is used to identify the layout
113 //____________________________________________________________________________
114 AliEMCAL::~AliEMCAL()
119 //____________________________________________________________________________
120 void AliEMCAL::InitConstants()
122 //initialize EMCAL values
124 fBirkC1 = 0.013/1.032;
125 fBirkC2 = 9.6e-6/(1.032 * 1.032);
127 fHighCharge = 8.2 ; // adjusted for a high gain range of 5.12 GeV (10 bits)
129 fHighLowGainFactor = 16. ; // adjusted for a low gain range of 82 GeV (10 bits)
130 fLowGainOffset = 1 ; // offset added to the module id to distinguish high and low gain data
133 //____________________________________________________________________________
134 void AliEMCAL::DefineMediumParameters()
137 // EMCAL cuts (Geant3)
139 Int_t * idtmed = fIdtmed->GetArray() - 1599 ;
140 // --- Set decent energy thresholds for gamma and electron tracking
142 // Tracking threshold for photons and electrons in Lead
143 Float_t cutgam=10.e-5; // 100 kev;
144 Float_t cutele=10.e-5; // 100 kev;
145 TString ntmp(GetTitle());
147 if(ntmp.Contains("10KEV")) {
148 cutele = cutgam = 1.e-5;
149 } else if(ntmp.Contains("50KEV")) {
150 cutele = cutgam = 5.e-5;
151 } else if(ntmp.Contains("100KEV")) {
152 cutele = cutgam = 1.e-4;
153 } else if(ntmp.Contains("200KEV")) {
154 cutele = cutgam = 2.e-4;
155 } else if(ntmp.Contains("500KEV")) {
156 cutele = cutgam = 5.e-4;
159 gMC->Gstpar(idtmed[1600],"CUTGAM", cutgam);
160 gMC->Gstpar(idtmed[1600],"CUTELE", cutele); // 1MEV -> 0.1MEV; 15-aug-05
161 gMC->Gstpar(idtmed[1600],"BCUTE", cutgam); // BCUTE and BCUTM start from GUTGUM
162 gMC->Gstpar(idtmed[1600],"BCUTM", cutgam); // BCUTE and BCUTM start from GUTGUM
163 // --- Generate explicitly delta rays in Lead ---
164 gMC->Gstpar(idtmed[1600], "LOSS", 3) ;
165 gMC->Gstpar(idtmed[1600], "DRAY", 1) ;
166 gMC->Gstpar(idtmed[1600], "DCUTE", cutele) ;
167 gMC->Gstpar(idtmed[1600], "DCUTM", cutele) ;
169 // --- in aluminium parts ---
170 gMC->Gstpar(idtmed[1602],"CUTGAM", cutgam) ;
171 gMC->Gstpar(idtmed[1602],"CUTELE", cutele) ;
172 gMC->Gstpar(idtmed[1602],"BCUTE", cutgam); // BCUTE and BCUTM start from GUTGUM
173 gMC->Gstpar(idtmed[1602],"BCUTM", cutgam); // BCUTE and BCUTM start from GUTGUM
174 gMC->Gstpar(idtmed[1602], "LOSS",3.) ;
175 gMC->Gstpar(idtmed[1602], "DRAY",1.) ;
176 gMC->Gstpar(idtmed[1602], "DCUTE", cutele) ;
177 gMC->Gstpar(idtmed[1602], "DCUTM", cutele) ;
179 // --- and finally thresholds for photons and electrons in the scintillator ---
180 gMC->Gstpar(idtmed[1601],"CUTGAM", cutgam) ;
181 gMC->Gstpar(idtmed[1601],"CUTELE", cutele) ;// 1MEV -> 0.1MEV; 15-aug-05
182 gMC->Gstpar(idtmed[1601],"BCUTE", cutgam); // BCUTE and BCUTM start from GUTGUM
183 gMC->Gstpar(idtmed[1601],"BCUTM", cutgam); // BCUTE and BCUTM start from GUTGUM
184 gMC->Gstpar(idtmed[1601], "LOSS",3) ; // generate delta rays
185 gMC->Gstpar(idtmed[1601], "DRAY",1) ;
186 gMC->Gstpar(idtmed[1601], "DCUTE", cutele) ;
187 gMC->Gstpar(idtmed[1601], "DCUTM", cutele) ;
190 gMC->Gstpar(idtmed[1603],"CUTGAM", cutgam);
191 gMC->Gstpar(idtmed[1603],"CUTELE", cutele);
192 gMC->Gstpar(idtmed[1603],"BCUTE", cutgam); // BCUTE and BCUTM start from GUTGUM
193 gMC->Gstpar(idtmed[1603],"BCUTM", cutgam); // BCUTE and BCUTM start from GUTGUM
194 // --- Generate explicitly delta rays
195 gMC->Gstpar(idtmed[1603], "LOSS",3);
196 gMC->Gstpar(idtmed[1603], "DRAY",1);
197 gMC->Gstpar(idtmed[1603], "DCUTE", cutele) ;
198 gMC->Gstpar(idtmed[1603], "DCUTM", cutele) ;
200 AliEMCALGeometry* geom = GetGeometry();
201 if(geom->GetILOSS()>=0) {
202 for(int i=1600; i<=1603; i++) gMC->Gstpar(idtmed[i], "LOSS", geom->GetILOSS()) ;
204 if(geom->GetIHADR()>=0) {
205 for(int i=1600; i<=1603; i++) gMC->Gstpar(idtmed[i], "HADR", geom->GetIHADR()) ;
209 //____________________________________________________________________________
210 AliDigitizer* AliEMCAL::CreateDigitizer(AliRunDigitizer* manager) const
212 //create and return the digitizer
213 return new AliEMCALDigitizer(manager);
216 //____________________________________________________________________________
217 void AliEMCAL::CreateMaterials()
219 // Definitions of materials to build EMCAL and associated tracking media.
220 // media number in idtmed are 1599 to 1698.
222 Float_t aAir[4]={12.0107,14.0067,15.9994,39.948};
223 Float_t zAir[4]={6.,7.,8.,18.};
224 Float_t wAir[4]={0.000124,0.755267,0.231781,0.012827};
225 Float_t dAir = 1.20479E-3;
226 AliMixture(0, "Air$", aAir, zAir, dAir, 4, wAir) ;
229 AliMaterial(1, "Pb$", 207.2, 82, 11.35, 0.56, 0., 0, 0) ;
232 // --- The polysterene scintillator (CH) ---
233 Float_t aP[2] = {12.011, 1.00794} ;
234 Float_t zP[2] = {6.0, 1.0} ;
235 Float_t wP[2] = {1.0, 1.0} ;
238 AliMixture(2, "Polystyrene$", aP, zP, dP, -2, wP) ;
241 AliMaterial(3, "Al$", 26.98, 13., 2.7, 8.9, 999., 0, 0) ;
242 // --- Absorption length is ignored ^
244 // 25-aug-04 by PAI - see PMD/AliPMDv0.cxx for STEEL definition
245 Float_t asteel[4] = { 55.847,51.9961,58.6934,28.0855 };
246 Float_t zsteel[4] = { 26.,24.,28.,14. };
247 Float_t wsteel[4] = { .715,.18,.1,.005 };
248 AliMixture(4, "STAINLESS STEEL$", asteel, zsteel, 7.88, 4, wsteel);
250 // DEFINITION OF THE TRACKING MEDIA
252 // for EMCAL: idtmed[1599->1698] equivalent to fIdtmed[0->100]
253 Int_t isxfld = gAlice->Field()->Integ() ;
254 Float_t sxmgmx = gAlice->Field()->Max() ;
256 // Air -> idtmed[1599]
257 AliMedium(0, "Air$", 0, 0,
258 isxfld, sxmgmx, 10.0, 1.0, 0.1, 0.1, 10.0, 0, 0) ;
260 // The Lead -> idtmed[1600]
262 AliMedium(1, "Lead$", 1, 0,
263 isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
265 // The scintillator of the CPV made of Polystyrene scintillator -> idtmed[1601]
266 float deemax = 0.1; // maximum fractional energy loss in one step (0 < DEEMAX â?¤ 1);i
267 AliMedium(2, "Scintillator$", 2, 1,
268 isxfld, sxmgmx, 10.0, 0.001, deemax, 0.001, 0.001, 0, 0) ;
270 // Various Aluminium parts made of Al -> idtmed[1602]
271 AliMedium(3, "Al$", 3, 0,
272 isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0) ;
274 // 25-aug-04 by PAI : see PMD/AliPMDv0.cxx for STEEL definition -> idtmed[1603]
275 AliMedium(4, "S steel$", 4, 0,
276 isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0) ;
279 //set constants for Birk's Law implentation
282 fBirkC2 = 9.6e-6/(dP * dP);
284 // Call just in case of Geant3; What to do in case of Geant4 ?
285 if(gMC->InheritsFrom("TGeant3")) DefineMediumParameters(); // Feb 20, 2007
287 //____________________________________________________________________________
288 void AliEMCAL::Digits2Raw()
290 // convert digits of the current event to raw data
292 AliEMCALLoader * loader = dynamic_cast<AliEMCALLoader*>(fLoader) ;
295 loader->LoadDigits("EMCAL");
297 TClonesArray* digits = loader->Digits() ;
300 Error("Digits2Raw", "no digits found !");
305 loader->LoadDigitizer();
306 AliEMCALDigitizer * digitizer = dynamic_cast<AliEMCALDigitizer *>(loader->Digitizer()) ;
309 AliEMCALGeometry* geom = GetGeometry();
311 AliError(Form("No geometry found !"));
315 AliAltroBuffer* buffer = NULL;
317 Int_t adcValuesLow[fgkTimeBins];
318 Int_t adcValuesHigh[fgkTimeBins];
320 //Load Mapping RCU files once
321 TString path = gSystem->Getenv("ALICE_ROOT");
322 path += "/EMCAL/mapping/RCU";
323 TString path0 = path+"0.data";//This file will change in future
324 TString path1 = path+"1.data";//This file will change in future
325 AliAltroMapping * mapping[2] ; // For the moment only 2
326 mapping[0] = new AliCaloAltroMapping(path0.Data());
327 mapping[1] = new AliCaloAltroMapping(path1.Data());
329 // loop over digits (assume ordered digits)
330 for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
331 AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
332 if (digit->GetAmp() < fgThreshold)
342 geom->GetCellIndex(digit->GetId(), nSM, nModule, nIphi, nIeta);
343 geom->GetCellPhiEtaIndexInSModule(nSM, nModule, nIphi, nIeta,iphi, ieta) ;
345 //Check which is the RCU of the cell.
348 if (0<=iphi&&iphi<8) iRCU=0; // first cable row
349 else if (8<=iphi&&iphi<16 && 0<=ieta&&ieta<24) iRCU=0; // first half;
352 else if(8<=iphi&&iphi<16 && 24<=ieta&&ieta<48) iRCU=1; // second half;
354 else if(16<=iphi&&iphi<24) iRCU=1; // third cable row
357 Int_t iDDL = fgDDLPerSuperModule* nSM + iRCU;
360 if (iDDL != prevDDL) {
361 // write real header and close previous file
365 buffer->WriteDataHeader(kFALSE, kFALSE);
369 // open new file and write dummy header
370 TString fileName = AliDAQ::DdlFileName("EMCAL",iDDL);
371 buffer = new AliAltroBuffer(fileName.Data(),mapping[iRCU]);
372 buffer->WriteDataHeader(kTRUE, kFALSE); //Dummy;
376 // out of time range signal (?)
377 if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
378 AliInfo("Signal is out of time range.\n");
379 buffer->FillBuffer((Int_t)digit->GetAmp());
380 buffer->FillBuffer(GetRawFormatTimeBins() ); // time bin
381 buffer->FillBuffer(3); // bunch length
382 buffer->WriteTrailer(3, ieta, iphi, nSM); // trailer
383 // calculate the time response function
386 Double_t energy = digit->GetAmp() * digitizer->GetECAchannel() + digitizer->GetECApedestal() ;
388 Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), energy, adcValuesHigh, adcValuesLow) ;
391 buffer->WriteChannel(ieta, iphi, 0, GetRawFormatTimeBins(), adcValuesLow, fgThreshold);
393 buffer->WriteChannel(ieta,iphi, 1, GetRawFormatTimeBins(), adcValuesHigh, fgThreshold);
398 // write real header and close last file
401 buffer->WriteDataHeader(kFALSE, kFALSE);
404 mapping[0]->Delete();
405 mapping[1]->Delete();
406 loader->UnloadDigits();
409 //____________________________________________________________________________
410 void AliEMCAL::Raw2Digits(AliRawReader* reader)
412 // convert raw data of the current event to digits
413 AliEMCALGeometry * geom = GetGeometry();
414 AliEMCALLoader * loader = dynamic_cast<AliEMCALLoader*>(fLoader) ;
417 loader->CleanDigits(); // start from scratch
418 loader->LoadDigits("EMCAL");
419 TClonesArray* digits = loader->Digits() ;
420 digits->Clear(); // yes, this is perhaps somewhat paranoid.. [clearing an extra time]
423 Error("Raw2Digits", "no digits found !");
427 Error("Raw2Digits", "no raw reader found !");
431 // and get the digitizer too
432 loader->LoadDigitizer();
433 AliEMCALDigitizer * digitizer = dynamic_cast<AliEMCALDigitizer *>(loader->Digitizer()) ;
435 // Use AliAltroRawStream to read the ALTRO format. No need to
436 // reinvent the wheel :-)
437 AliCaloRawStream in(reader,"EMCAL");
438 // Select EMCAL DDL's;
439 reader->Select("EMCAL");
441 // reading is from previously existing AliEMCALGetter.cxx
443 Bool_t first = kTRUE ;
445 TF1 * signalF = new TF1("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
446 signalF->SetParNames("Charge", "Gain", "Amplitude", "TimeZero");
448 Bool_t lowGainFlag = kFALSE ;
453 Double_t energy = 0. ;
455 TGraph * gLowGain = new TGraph(GetRawFormatTimeBins()) ;
456 TGraph * gHighGain= new TGraph(GetRawFormatTimeBins()) ;
458 while ( in.Next() ) { // EMCAL entries loop
459 if ( in.IsNewRow() ) {//phi
460 if ( in.IsNewColumn() ) {//eta
461 id = geom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
463 FitRaw(lowGainFlag, gLowGain, gHighGain, signalF, energy, time) ;
465 if (time == 0. && energy == 0.) {
469 amp = static_cast<Int_t>( (energy - digitizer->GetECApedestal()) / digitizer->GetECAchannel() + 0.5 ) ;
473 new((*digits)[idigit]) AliEMCALDigit( -1, -1, id, amp, time) ;
477 for (index = 0; index < GetRawFormatTimeBins(); index++) {
478 gLowGain->SetPoint(index, index * GetRawFormatTimeMax() / GetRawFormatTimeBins(), 0) ;
479 gHighGain->SetPoint(index, index * GetRawFormatTimeMax() / GetRawFormatTimeBins(), 0) ;
483 id = geom->GetAbsCellIdFromCellIndexes(in.GetModule(), in.GetRow(), in.GetColumn()) ;
484 if (in.GetModule() == GetRawFormatLowGainOffset() ) {
485 lowGainFlag = kTRUE ;
488 lowGainFlag = kFALSE ;
493 gLowGain->SetPoint(in.GetTime(),
494 in.GetTime()* GetRawFormatTimeMax() / GetRawFormatTimeBins(),
498 gHighGain->SetPoint(in.GetTime(),
499 in.GetTime() * GetRawFormatTimeMax() / GetRawFormatTimeBins(),
502 } // EMCAL entries loop
512 //____________________________________________________________________________
513 void AliEMCAL::FitRaw(Bool_t lowGainFlag, TGraph * gLowGain, TGraph * gHighGain, TF1* signalF, Double_t & energy, Double_t & time)
515 // Fits the raw signal time distribution; from AliEMCALGetter
517 const Int_t kNoiseThreshold = 0 ;
518 Double_t timezero1 = 0., timezero2 = 0., timemax = 0. ;
519 Double_t signal = 0., signalmax = 0. ;
523 timezero1 = timezero2 = signalmax = timemax = 0. ;
524 signalF->FixParameter(0, GetRawFormatLowCharge()) ;
525 signalF->FixParameter(1, GetRawFormatLowGain()) ;
527 for (index = 0; index < GetRawFormatTimeBins(); index++) {
528 gLowGain->GetPoint(index, time, signal) ;
529 if (signal > kNoiseThreshold && timezero1 == 0.)
531 if (signal <= kNoiseThreshold && timezero1 > 0. && timezero2 == 0.)
533 if (signal > signalmax) {
538 signalmax /= RawResponseFunctionMax(GetRawFormatLowCharge(),
539 GetRawFormatLowGain()) ;
540 if ( timezero1 + GetRawFormatTimePeak() < GetRawFormatTimeMax() * 0.4 ) { // else its noise
541 signalF->SetParameter(2, signalmax) ;
542 signalF->SetParameter(3, timezero1) ;
543 gLowGain->Fit(signalF, "QRON", "", 0., timezero2); //, "QRON") ;
544 energy = signalF->GetParameter(2) ;
545 time = signalF->GetMaximumX() - GetRawFormatTimePeak() - GetRawFormatTimeTrigger() ;
548 timezero1 = timezero2 = signalmax = timemax = 0. ;
549 signalF->FixParameter(0, GetRawFormatHighCharge()) ;
550 signalF->FixParameter(1, GetRawFormatHighGain()) ;
552 for (index = 0; index < GetRawFormatTimeBins(); index++) {
553 gHighGain->GetPoint(index, time, signal) ;
554 if (signal > kNoiseThreshold && timezero1 == 0.)
556 if (signal <= kNoiseThreshold && timezero1 > 0. && timezero2 == 0.)
558 if (signal > signalmax) {
563 signalmax /= RawResponseFunctionMax(GetRawFormatHighCharge(),
564 GetRawFormatHighGain()) ;;
565 if ( timezero1 + GetRawFormatTimePeak() < GetRawFormatTimeMax() * 0.4 ) { // else its noise
566 signalF->SetParameter(2, signalmax) ;
567 signalF->SetParameter(3, timezero1) ;
568 gHighGain->Fit(signalF, "QRON", "", 0., timezero2) ;
569 energy = signalF->GetParameter(2) ;
570 time = signalF->GetMaximumX() - GetRawFormatTimePeak() - GetRawFormatTimeTrigger() ;
577 //____________________________________________________________________________
578 void AliEMCAL::Hits2SDigits()
580 // create summable digits
583 AliEMCALSDigitizer emcalDigitizer(fLoader->GetRunLoader()->GetFileName().Data()) ;
584 emcalDigitizer.SetEventRange(0, -1) ; // do all the events
585 emcalDigitizer.ExecuteTask() ;
588 //____________________________________________________________________________
590 AliLoader* AliEMCAL::MakeLoader(const char* topfoldername)
592 //different behaviour than standard (singleton getter)
593 // --> to be discussed and made eventually coherent
594 fLoader = new AliEMCALLoader(GetName(),topfoldername);
598 //__________________________________________________________________
599 Double_t AliEMCAL::RawResponseFunction(Double_t *x, Double_t *par)
601 // Shape of the electronics raw reponse:
602 // It is a semi-gaussian, 2nd order Gamma function of the general form
603 // v(t) = n**n * Q * A**n / C *(t/tp)**n * exp(-n * t/tp) with
604 // tp : peaking time par[0]
605 // n : order of the function
606 // C : integrating capacitor in the preamplifier
607 // A : open loop gain of the preamplifier
608 // Q : the total APD charge to be measured Q = C * energy
611 Double_t xx = x[0] - ( fgTimeTrigger + par[3] ) ;
613 if (xx < 0 || xx > fgTimeMax)
616 Double_t fac = par[0] * TMath::Power(fgOrder, fgOrder) * TMath::Power(par[1], fgOrder) / fgCapa ;
617 signal = fac * par[2] * TMath::Power(xx / fgTimePeak, fgOrder) * TMath::Exp(-fgOrder * (xx / fgTimePeak)) ;
622 //__________________________________________________________________
623 Double_t AliEMCAL::RawResponseFunctionMax(Double_t charge, Double_t gain)
625 //compute the maximum of the raw response function and return
626 return ( charge * TMath::Power(fgOrder, fgOrder) * TMath::Power(gain, fgOrder)
627 / ( fgCapa * TMath::Exp(fgOrder) ) );
630 //__________________________________________________________________
631 Bool_t AliEMCAL::RawSampledResponse(
632 const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL) const
634 // for a start time dtime and an amplitude damp given by digit,
635 // calculates the raw sampled response AliEMCAL::RawResponseFunction
637 const Int_t kRawSignalOverflow = 0x3FF ;
638 Bool_t lowGain = kFALSE ;
640 TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
642 for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
643 signalF.SetParameter(0, GetRawFormatHighCharge() ) ;
644 signalF.SetParameter(1, GetRawFormatHighGain() ) ;
645 signalF.SetParameter(2, damp) ;
646 signalF.SetParameter(3, dtime) ;
647 Double_t time = iTime * GetRawFormatTimeMax() / GetRawFormatTimeBins() ;
648 Double_t signal = signalF.Eval(time) ;
649 if ( static_cast<Int_t>(signal+0.5) > kRawSignalOverflow ){ // larger than 10 bits
650 signal = kRawSignalOverflow ;
653 adcH[iTime] = static_cast<Int_t>(signal + 0.5) ;
655 signalF.SetParameter(0, GetRawFormatLowCharge() ) ;
656 signalF.SetParameter(1, GetRawFormatLowGain() ) ;
657 signal = signalF.Eval(time) ;
658 if ( static_cast<Int_t>(signal+0.5) > kRawSignalOverflow) // larger than 10 bits
659 signal = kRawSignalOverflow ;
660 adcL[iTime] = static_cast<Int_t>(0.5 + signal ) ;