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