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