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