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