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