]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCAL.cxx
Major changes in AliAltroBuffer. Now it can be used only for writing of raw data...
[u/mrichter/AliRoot.git] / EMCAL / AliEMCAL.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
18 //_________________________________________________________________________
19 // Base Class for EMCAL description:
20 // This class contains material definitions    
21 // for the EMCAL - It does not place the detector in Alice
22 //*-- Author: Yves Schutz (SUBATECH) 
23 //
24 //*-- Additional Contributions: Sahal Yacoob (LBNL/UCT)
25 //
26 //////////////////////////////////////////////////////////////////////////////
27
28 // --- ROOT system ---
29 class TFile;
30 #include <TFolder.h> 
31 #include <TTree.h>
32 #include <TVirtualMC.h> 
33 #include <TH1F.h> 
34 #include <TF1.h> 
35 #include <TRandom.h> 
36 #include <TGraph.h> 
37
38 // --- Standard library ---
39
40 // --- AliRoot header files ---
41 #include "AliMagF.h"
42 #include "AliEMCAL.h"
43 #include "AliRun.h"
44 #include "AliEMCALLoader.h"
45 #include "AliEMCALSDigitizer.h"
46 #include "AliEMCALDigitizer.h"
47 #include "AliEMCALDigit.h"
48 #include "AliAltroBuffer.h"
49 #include "AliRawReader.h"
50 #include "AliEMCALRawStream.h"
51
52 ClassImp(AliEMCAL)
53 Double_t AliEMCAL::fgCapa        = 1.;        // 1pF 
54 Int_t    AliEMCAL::fgOrder       = 2 ;
55 Double_t AliEMCAL::fgTimeMax     = 2.56E-5 ;  // each sample is over 100 ns fTimeMax/fTimeBins
56 Double_t AliEMCAL::fgTimePeak    = 4.1E-6 ;   // 4 micro seconds
57 Double_t AliEMCAL::fgTimeTrigger = 100E-9 ;      // 100ns, just for a reference
58 // some digitization constants
59 Int_t    AliEMCAL::fgDDLOffset = 0x800;
60 Int_t    AliEMCAL::fgThreshold = 1;
61 // 24*48=1152 towers per SM; divided up on 3 DDLs, 
62 // each DDL with 12FEC *32towers or 12*32*2 channels (high&low gain) 
63 Int_t    AliEMCAL::fgChannelsPerDDL = 768; // 2*(1152/3 or 12*32) 
64  
65 //____________________________________________________________________________
66 AliEMCAL::AliEMCAL():AliDetector()
67 {
68   // Default ctor 
69   fName = "EMCAL" ;
70 }
71
72 //____________________________________________________________________________
73 AliEMCAL::AliEMCAL(const char* name, const char* title): AliDetector(name,title)
74 {
75   //   ctor : title is used to identify the layout
76
77   fHighCharge        = 8.2 ;          // adjusted for a high gain range of 5.12 GeV (10 bits)
78   fHighGain          = 6.64 ; 
79   fHighLowGainFactor = 16. ;          // adjusted for a low gain range of 82 GeV (10 bits) 
80   fLowGainOffset     = 1 ;            // offset added to the module id to distinguish high and low gain data
81 }
82
83 //____________________________________________________________________________
84 AliEMCAL::~AliEMCAL()
85 {
86   //dtor
87 }
88
89 //____________________________________________________________________________
90 void AliEMCAL::Copy(AliEMCAL & emcal) const
91 {
92   //copy
93
94   TObject::Copy(emcal) ; 
95   emcal.fHighCharge        = fHighCharge ;
96   emcal.fHighGain          = fHighGain ; 
97   emcal.fHighLowGainFactor = fHighLowGainFactor ;  
98   emcal.fLowGainOffset     = fLowGainOffset;   
99 }
100
101 //____________________________________________________________________________
102 AliDigitizer* AliEMCAL::CreateDigitizer(AliRunDigitizer* manager) const
103 {
104   //create and return the digitizer
105   return new AliEMCALDigitizer(manager);
106 }
107
108 //____________________________________________________________________________
109 void AliEMCAL::CreateMaterials()
110 {
111   // Definitions of materials to build EMCAL and associated tracking media.
112   // media number in idtmed are 1599 to 1698.
113
114   // --- Air ---               
115   Float_t aAir[4]={12.0107,14.0067,15.9994,39.948};
116   Float_t zAir[4]={6.,7.,8.,18.};
117   Float_t wAir[4]={0.000124,0.755267,0.231781,0.012827};
118   Float_t dAir = 1.20479E-3;
119   AliMixture(0, "Air$", aAir, zAir, dAir, 4, wAir) ;
120
121   // --- Lead ---                                                                     
122   AliMaterial(1, "Pb$", 207.2, 82, 11.35, 0.56, 0., 0, 0) ;
123
124
125   // --- The polysterene scintillator (CH) ---
126   Float_t aP[2] = {12.011, 1.00794} ;
127   Float_t zP[2] = {6.0, 1.0} ;
128   Float_t wP[2] = {1.0, 1.0} ;
129   Float_t dP = 1.032 ;
130
131   AliMixture(2, "Polystyrene$", aP, zP, dP, -2, wP) ;
132
133   // --- Aluminium ---
134   AliMaterial(3, "Al$", 26.98, 13., 2.7, 8.9, 999., 0, 0) ;
135   // ---         Absorption length is ignored ^
136
137   // 25-aug-04 by PAI - see  PMD/AliPMDv0.cxx for STEEL definition
138   Float_t asteel[4] = { 55.847,51.9961,58.6934,28.0855 };
139   Float_t zsteel[4] = { 26.,24.,28.,14. };
140   Float_t wsteel[4] = { .715,.18,.1,.005 };
141   AliMixture(4, "STAINLESS STEEL$", asteel, zsteel, 7.88, 4, wsteel);
142
143   // DEFINITION OF THE TRACKING MEDIA
144
145   // for EMCAL: idtmed[1599->1698] equivalent to fIdtmed[0->100]
146   Int_t * idtmed = fIdtmed->GetArray() - 1599 ; 
147   Int_t   isxfld = gAlice->Field()->Integ() ;
148   Float_t sxmgmx = gAlice->Field()->Max() ;
149
150   // Air                                                                         -> idtmed[1599]
151  AliMedium(0, "Air$", 0, 0,
152              isxfld, sxmgmx, 10.0, 1.0, 0.1, 0.1, 10.0, 0, 0) ;
153
154   // The Lead                                                                      -> idtmed[1600]
155  
156   AliMedium(1, "Lead$", 1, 0,
157              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0) ;
158
159  // The scintillator of the CPV made of Polystyrene scintillator                   -> idtmed[1601]
160   AliMedium(2, "Scintillator$", 2, 1,
161             isxfld, sxmgmx, 10.0, 0.001, 0.1, 0.001, 0.001, 0, 0) ;
162
163   // Various Aluminium parts made of Al                                            -> idtmed[1602]
164   AliMedium(3, "Al$", 3, 0,
165              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0) ;
166
167   // 25-aug-04 by PAI : see  PMD/AliPMDv0.cxx for STEEL definition                 -> idtmed[1603]
168   AliMedium(4, "S steel$", 4, 0, 
169              isxfld, sxmgmx, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0) ;
170
171 // --- Set decent energy thresholds for gamma and electron tracking
172
173   // Tracking threshold for photons and electrons in Lead 
174   Float_t cutgam=10.e-5; // 100 kev;
175   Float_t cutele=10.e-5; // 100 kev;
176   TString ntmp(GetTitle()); 
177   ntmp.ToUpper();
178   if(ntmp.Contains("10KEV")) {
179     cutele = cutgam = 1.e-5;
180   } else if(ntmp.Contains("50KEV")) {
181     cutele = cutgam = 5.e-5;
182   } else if(ntmp.Contains("100KEV")) {
183     cutele = cutgam = 1.e-4;
184   } else if(ntmp.Contains("200KEV")) {
185     cutele = cutgam = 2.e-4;
186   } else if(ntmp.Contains("500KEV")) {
187     cutele = cutgam = 5.e-4;
188   }
189
190   gMC->Gstpar(idtmed[1600],"CUTGAM", cutgam);
191   gMC->Gstpar(idtmed[1600],"CUTELE", cutele); // 1MEV -> 0.1MEV; 15-aug-05
192   gMC->Gstpar(idtmed[1600],"BCUTE",  cutgam);  // BCUTE and BCUTM start from GUTGUM
193   gMC->Gstpar(idtmed[1600],"BCUTM",  cutgam);  // BCUTE and BCUTM start from GUTGUM
194   // --- Generate explicitly delta rays in Lead ---
195   gMC->Gstpar(idtmed[1600], "LOSS",3.) ;
196   gMC->Gstpar(idtmed[1600], "DRAY",1.) ;
197   gMC->Gstpar(idtmed[1600], "DCUTE", cutele) ;
198   gMC->Gstpar(idtmed[1600], "DCUTM", cutele) ;
199
200 // --- in aluminium parts ---
201   gMC->Gstpar(idtmed[1602],"CUTGAM", cutgam) ;
202   gMC->Gstpar(idtmed[1602],"CUTELE", cutele) ;
203   gMC->Gstpar(idtmed[1602],"BCUTE",  cutgam);  // BCUTE and BCUTM start from GUTGUM
204   gMC->Gstpar(idtmed[1602],"BCUTM",  cutgam);  // BCUTE and BCUTM start from GUTGUM
205   gMC->Gstpar(idtmed[1602], "LOSS",3.) ;
206   gMC->Gstpar(idtmed[1602], "DRAY",1.) ;
207   gMC->Gstpar(idtmed[1602], "DCUTE", cutele) ;
208   gMC->Gstpar(idtmed[1602], "DCUTM", cutele) ;
209
210 // --- and finally thresholds for photons and electrons in the scintillator ---
211   gMC->Gstpar(idtmed[1601],"CUTGAM", cutgam) ;
212   gMC->Gstpar(idtmed[1601],"CUTELE", cutele) ;// 1MEV -> 0.1MEV; 15-aug-05
213   gMC->Gstpar(idtmed[1601],"BCUTE",  cutgam);  // BCUTE and BCUTM start from GUTGUM
214   gMC->Gstpar(idtmed[1601],"BCUTM",  cutgam);  // BCUTE and BCUTM start from GUTGUM
215   gMC->Gstpar(idtmed[1601], "LOSS",3.) ; // generate delta rays 
216   gMC->Gstpar(idtmed[1601], "DRAY",1.) ;
217   gMC->Gstpar(idtmed[1601], "DCUTE", cutele) ;
218   gMC->Gstpar(idtmed[1601], "DCUTM", cutele) ;
219
220   // S steel - 
221   gMC->Gstpar(idtmed[1603],"CUTGAM", cutgam);
222   gMC->Gstpar(idtmed[1603],"CUTELE", cutele);
223   gMC->Gstpar(idtmed[1603],"BCUTE",  cutgam);  // BCUTE and BCUTM start from GUTGUM
224   gMC->Gstpar(idtmed[1603],"BCUTM",  cutgam);  // BCUTE and BCUTM start from GUTGUM
225   // --- Generate explicitly delta rays 
226   gMC->Gstpar(idtmed[1603], "LOSS",3.);
227   gMC->Gstpar(idtmed[1603], "DRAY",1.);
228   gMC->Gstpar(idtmed[1603], "DCUTE", cutele) ;
229   gMC->Gstpar(idtmed[1603], "DCUTM", cutele) ;
230
231   //set constants for Birk's Law implentation
232   fBirkC0 =  1;
233   fBirkC1 =  0.013/dP;
234   fBirkC2 =  9.6e-6/(dP * dP);
235
236 }
237       
238 //____________________________________________________________________________
239 void AliEMCAL::Digits2Raw()
240 {
241   // convert digits of the current event to raw data
242   AliEMCALLoader * loader = dynamic_cast<AliEMCALLoader*>(fLoader) ; 
243
244   // get the digits
245   loader->LoadDigits("EMCAL");
246   TClonesArray* digits = loader->Digits() ;
247
248   if (!digits) {
249     Error("Digits2Raw", "no digits found !");
250     return;
251   }
252
253   // get the digitizer 
254   loader->LoadDigitizer();
255   AliEMCALDigitizer * digitizer = dynamic_cast<AliEMCALDigitizer *>(loader->Digitizer())  ; 
256   
257
258   AliAltroBuffer* buffer = NULL;
259   Int_t prevDDL = -1;
260   Int_t adcValuesLow[fgkTimeBins];
261   Int_t adcValuesHigh[fgkTimeBins];
262   
263   // loop over digits (assume ordered digits)
264   for (Int_t iDigit = 0; iDigit < digits->GetEntries(); iDigit++) {
265     AliEMCALDigit* digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit)) ;
266     if (digit->GetAmp() < fgThreshold) 
267       continue;
268     Int_t iDDL = digit->GetId() / fgChannelsPerDDL ;
269     // for each DDL id is numbered from 1 to  fgChannelsperDDL -1 
270     Int_t idDDL = digit->GetId() - iDDL * ( fgChannelsPerDDL - 1 ) ;  
271     // new DDL
272     if (iDDL != prevDDL) {
273       // write real header and close previous file
274       if (buffer) {
275         buffer->Flush();
276         buffer->WriteDataHeader(kFALSE, kFALSE);
277         delete buffer;
278       }
279
280       // open new file and write dummy header
281       TString fileName("EMCAL_") ;
282       fileName += (iDDL + fgDDLOffset) ; 
283       fileName += ".ddl" ; 
284       buffer = new AliAltroBuffer(fileName.Data());
285       buffer->WriteDataHeader(kTRUE, kFALSE);  //Dummy;
286
287       prevDDL = iDDL;
288     }
289
290     // out of time range signal (?)
291     if (digit->GetTimeR() > GetRawFormatTimeMax() ) {
292       buffer->FillBuffer(digit->GetAmp());
293       buffer->FillBuffer(GetRawFormatTimeBins() );  // time bin
294       buffer->FillBuffer(3);          // bunch length
295       buffer->WriteTrailer(3, idDDL, 0, 0);  // trailer
296
297       // calculate the time response function
298     } else {
299       Double_t energy = 0 ;  
300       energy = digit->GetAmp() * digitizer->GetECAchannel() + digitizer->GetECApedestal() ; 
301       
302       Bool_t lowgain = RawSampledResponse(digit->GetTimeR(), energy, adcValuesHigh, adcValuesLow) ; 
303       
304       if (lowgain) 
305         buffer->WriteChannel(iDDL, 0, fLowGainOffset, 
306                              GetRawFormatTimeBins(), adcValuesLow, fgThreshold);
307       else 
308         buffer->WriteChannel(iDDL, 0, 0, 
309                              GetRawFormatTimeBins(), adcValuesHigh, fgThreshold);
310       
311     }
312   }
313   
314   // write real header and close last file
315   if (buffer) {
316     buffer->Flush();
317     buffer->WriteDataHeader(kFALSE, kFALSE);
318     delete buffer;
319   }
320
321   loader->UnloadDigits();
322 }
323
324 //____________________________________________________________________________
325 void AliEMCAL::Raw2Digits(AliRawReader* reader)
326 {
327   // convert raw data of the current event to digits
328   GetGeometry();
329   AliEMCALLoader * loader = dynamic_cast<AliEMCALLoader*>(fLoader) ; 
330
331   // get the digits
332   loader->CleanDigits(); // start from scratch
333   loader->LoadDigits("EMCAL");
334   TClonesArray* digits = loader->Digits() ;
335   digits->Clear(); // yes, this is perhaps somewhat paranoid.. [clearing an extra time]
336
337   if (!digits) {
338     Error("Raw2Digits", "no digits found !");
339     return;
340   }
341   if (!reader) {
342     Error("Raw2Digits", "no raw reader found !");
343     return;
344   }
345
346   // and get the digitizer too 
347   loader->LoadDigitizer();
348   AliEMCALDigitizer * digitizer = dynamic_cast<AliEMCALDigitizer *>(loader->Digitizer())  ; 
349
350   // Use AliAltroRawStream to read the ALTRO format.  No need to
351   // reinvent the wheel :-) 
352   AliEMCALRawStream in(reader);
353   // Select EMCAL DDL's; lowest 8 bits of DDL offser is used for something else.. 
354   reader->Select(fgDDLOffset >> 8);
355
356   // reading is from previously existing AliEMCALGetter.cxx
357   // ReadRaw method
358   Bool_t first = kTRUE ;
359  
360   TF1 * signalF = new TF1("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
361   signalF->SetParNames("Charge", "Gain", "Amplitude", "TimeZero"); 
362   
363   Int_t id = -1;
364   Bool_t lowGainFlag = kFALSE ; 
365
366   Int_t idigit = 0 ; 
367   Int_t amp = 0 ; 
368   Double_t time = 0. ; 
369   Double_t energy = 0. ; 
370
371   TGraph * gLowGain = new TGraph(GetRawFormatTimeBins()) ; 
372   TGraph * gHighGain= new TGraph(GetRawFormatTimeBins()) ;  
373
374   while ( in.Next() ) { // EMCAL entries loop 
375     if ( in.IsNewId() ) {
376       if (!first) {
377         FitRaw(lowGainFlag, gLowGain, gHighGain, signalF, energy, time) ; 
378
379         if (time == 0. && energy == 0.) { 
380           amp = 0 ; 
381         }
382         else {
383           amp = static_cast<Int_t>( (energy - digitizer->GetECApedestal()) / digitizer->GetECAchannel() + 0.5 ) ; 
384         }
385
386         if (amp > 0) {
387           new((*digits)[idigit]) AliEMCALDigit( -1, -1, id, amp, time) ;        
388           idigit++ ; 
389         }
390         Int_t index ; 
391         for (index = 0; index < GetRawFormatTimeBins(); index++) {
392           gLowGain->SetPoint(index, index * GetRawFormatTimeMax() / GetRawFormatTimeBins(), 0) ;  
393           gHighGain->SetPoint(index, index * GetRawFormatTimeMax() / GetRawFormatTimeBins(), 0) ; 
394         } 
395       } // not first  
396       first = kFALSE ; 
397       id = in.GetId() ; 
398       if (in.GetModule() == GetRawFormatLowGainOffset() ) {
399         lowGainFlag = kTRUE ; 
400       }
401       else { 
402         lowGainFlag = kFALSE ; 
403       }
404     } // new Id?
405     if (lowGainFlag) {
406       gLowGain->SetPoint(in.GetTime(), 
407                          in.GetTime()* GetRawFormatTimeMax() / GetRawFormatTimeBins(), 
408                          in.GetSignal()) ;
409     }
410     else { 
411       gHighGain->SetPoint(in.GetTime(), 
412                           in.GetTime() * GetRawFormatTimeMax() / GetRawFormatTimeBins(), 
413                           in.GetSignal() ) ;
414     }
415   } // EMCAL entries loop
416   digits->Sort() ; 
417
418   delete signalF ; 
419   delete gLowGain;
420   delete gHighGain ; 
421     
422   return ; 
423 }
424
425 //____________________________________________________________________________ 
426 void AliEMCAL::FitRaw(Bool_t lowGainFlag, TGraph * gLowGain, TGraph * gHighGain, TF1* signalF, Double_t & energy, Double_t & time)
427 {
428   // Fits the raw signal time distribution; from AliEMCALGetter 
429
430   const Int_t kNoiseThreshold = 0 ;
431   Double_t timezero1 = 0., timezero2 = 0., timemax = 0. ;
432   Double_t signal = 0., signalmax = 0. ;       
433   energy = time = 0. ; 
434
435   if (lowGainFlag) {
436     timezero1 = timezero2 = signalmax = timemax = 0. ;
437     signalF->FixParameter(0, GetRawFormatLowCharge()) ; 
438     signalF->FixParameter(1, GetRawFormatLowGain()) ; 
439     Int_t index ; 
440     for (index = 0; index < GetRawFormatTimeBins(); index++) {
441       gLowGain->GetPoint(index, time, signal) ; 
442       if (signal > kNoiseThreshold && timezero1 == 0.) 
443         timezero1 = time ;
444       if (signal <= kNoiseThreshold && timezero1 > 0. && timezero2 == 0.)
445         timezero2 = time ; 
446       if (signal > signalmax) {
447         signalmax = signal ; 
448         timemax   = time ; 
449       }
450     }
451     signalmax /= RawResponseFunctionMax(GetRawFormatLowCharge(), 
452                                                 GetRawFormatLowGain()) ;
453     if ( timezero1 + GetRawFormatTimePeak() < GetRawFormatTimeMax() * 0.4 ) { // else its noise 
454       signalF->SetParameter(2, signalmax) ; 
455       signalF->SetParameter(3, timezero1) ;         
456       gLowGain->Fit(signalF, "QRON", "", 0., timezero2); //, "QRON") ; 
457       energy = signalF->GetParameter(2) ; 
458       time   = signalF->GetMaximumX() - GetRawFormatTimePeak() - GetRawFormatTimeTrigger() ;
459     }
460   } else {
461     timezero1 = timezero2 = signalmax = timemax = 0. ;
462     signalF->FixParameter(0, GetRawFormatHighCharge()) ; 
463     signalF->FixParameter(1, GetRawFormatHighGain()) ; 
464     Int_t index ; 
465     for (index = 0; index < GetRawFormatTimeBins(); index++) {
466       gHighGain->GetPoint(index, time, signal) ;               
467       if (signal > kNoiseThreshold && timezero1 == 0.) 
468         timezero1 = time ;
469       if (signal <= kNoiseThreshold && timezero1 > 0. && timezero2 == 0.)
470         timezero2 = time ; 
471       if (signal > signalmax) {
472         signalmax = signal ;   
473         timemax   = time ; 
474       }
475     }
476     signalmax /= RawResponseFunctionMax(GetRawFormatHighCharge(), 
477                                                 GetRawFormatHighGain()) ;;
478     if ( timezero1 + GetRawFormatTimePeak() < GetRawFormatTimeMax() * 0.4 ) { // else its noise  
479       signalF->SetParameter(2, signalmax) ; 
480       signalF->SetParameter(3, timezero1) ;               
481       gHighGain->Fit(signalF, "QRON", "", 0., timezero2) ; 
482       energy = signalF->GetParameter(2) ; 
483       time   = signalF->GetMaximumX() - GetRawFormatTimePeak() - GetRawFormatTimeTrigger() ;
484     }
485   }
486   
487   return;
488 }
489
490 //____________________________________________________________________________
491 void AliEMCAL::Hits2SDigits()  
492
493 // create summable digits
494
495   GetGeometry();
496   AliEMCALSDigitizer emcalDigitizer(fLoader->GetRunLoader()->GetFileName().Data()) ;
497   emcalDigitizer.SetEventRange(0, -1) ; // do all the events
498   emcalDigitizer.ExecuteTask() ;
499 }
500
501 //____________________________________________________________________________
502
503 AliLoader* AliEMCAL::MakeLoader(const char* topfoldername)
504 {
505 //different behaviour than standard (singleton getter)
506 // --> to be discussed and made eventually coherent
507  fLoader = new AliEMCALLoader(GetName(),topfoldername);
508  return fLoader;
509 }
510
511 //__________________________________________________________________
512 Double_t AliEMCAL::RawResponseFunction(Double_t *x, Double_t *par)
513 {
514   // Shape of the electronics raw reponse:
515   // It is a semi-gaussian, 2nd order Gamma function of the general form
516   // v(t) = n**n * Q * A**n / C *(t/tp)**n * exp(-n * t/tp) with 
517   // tp : peaking time par[0]
518   // n  : order of the function
519   // C  : integrating capacitor in the preamplifier
520   // A  : open loop gain of the preamplifier
521   // Q  : the total APD charge to be measured Q = C * energy
522   
523   Double_t signal ;
524   Double_t xx = x[0] - ( fgTimeTrigger + par[3] ) ; 
525   
526   if (xx < 0 || xx > fgTimeMax) 
527     signal = 0. ;  
528   else { 
529     Double_t fac = par[0] * TMath::Power(fgOrder, fgOrder) * TMath::Power(par[1], fgOrder) / fgCapa ; 
530     signal = fac * par[2] * TMath::Power(xx / fgTimePeak, fgOrder) * TMath::Exp(-fgOrder * (xx / fgTimePeak)) ; 
531   }
532   return signal ;  
533 }
534
535 //__________________________________________________________________
536 Double_t AliEMCAL::RawResponseFunctionMax(Double_t charge, Double_t gain) 
537 {
538   //compute the maximum of the raw response function and return
539   return ( charge * TMath::Power(fgOrder, fgOrder) * TMath::Power(gain, fgOrder) 
540      / ( fgCapa * TMath::Exp(fgOrder) ) );  
541
542 }
543 //__________________________________________________________________
544 Bool_t AliEMCAL::RawSampledResponse(
545 const Double_t dtime, const Double_t damp, Int_t * adcH, Int_t * adcL) const 
546 {
547   // for a start time dtime and an amplitude damp given by digit, 
548   // calculates the raw sampled response AliEMCAL::RawResponseFunction
549
550   const Int_t kRawSignalOverflow = 0x3FF ; 
551   Bool_t lowGain = kFALSE ; 
552
553   TF1 signalF("signal", RawResponseFunction, 0, GetRawFormatTimeMax(), 4);
554
555   for (Int_t iTime = 0; iTime < GetRawFormatTimeBins(); iTime++) {
556     signalF.SetParameter(0, GetRawFormatHighCharge() ) ; 
557     signalF.SetParameter(1, GetRawFormatHighGain() ) ; 
558     signalF.SetParameter(2, damp) ; 
559     signalF.SetParameter(3, dtime) ; 
560     Double_t time = iTime * GetRawFormatTimeMax() / GetRawFormatTimeBins() ;
561     Double_t signal = signalF.Eval(time) ;     
562     if ( static_cast<Int_t>(signal+0.5) > kRawSignalOverflow ){  // larger than 10 bits 
563       signal = kRawSignalOverflow ;
564       lowGain = kTRUE ; 
565     }
566     adcH[iTime] =  static_cast<Int_t>(signal + 0.5) ;
567
568     signalF.SetParameter(0, GetRawFormatLowCharge() ) ;     
569     signalF.SetParameter(1, GetRawFormatLowGain() ) ; 
570     signal = signalF.Eval(time) ;  
571     if ( static_cast<Int_t>(signal+0.5) > kRawSignalOverflow)  // larger than 10 bits 
572       signal = kRawSignalOverflow ;
573     adcL[iTime] = static_cast<Int_t>(0.5 + signal ) ; 
574
575   }
576   return lowGain ; 
577 }