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