]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALDigitizer.cxx
Changing to centrality dependent corrections
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALDigitizer.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
18 //_________________________________________________________________________
19 // 
20 //////////////////////////////////////////////////////////////////////////////
21 // Class performs digitization of Summable digits from simulated data
22 //  
23 // In addition it performs mixing/embedding of summable digits from different events.
24 //
25 // For each event 3 branches are created in TreeD:
26 //   "EMCAL" - list of digits
27 //   "EMCALTRG" - list of trigger digits
28 //   "AliEMCALDigitizer" - AliEMCALDigitizer with all parameters used in digitization
29 //
30 //
31 ////////////////////////////////////////////////////////////////////////////////////
32 //
33 //*-- Author: Sahal Yacoob (LBL)
34 // based on : AliEMCALDigitizer
35 // Modif: 
36 //  August 2002 Yves Schutz: clone PHOS as closely as possible and intoduction
37 //                           of new IO (a la PHOS)
38 //  November 2003 Aleksei Pavlinov : adopted for Shish-Kebab geometry 
39 //  July 2011 GCB: Digitizer modified to accomodate embedding. 
40 //                 Time calibration added. Decalibration possibility of energy and time added
41 //_________________________________________________________________________________
42
43 // --- ROOT system ---
44 #include <TROOT.h>
45 #include <TTree.h>
46 #include <TSystem.h>
47 #include <TBenchmark.h>
48 #include <TBrowser.h>
49 #include <TObjectTable.h>
50 #include <TRandom.h>
51 #include <TF1.h>
52 #include <cassert>
53
54 // --- AliRoot header files ---
55 #include "AliLog.h"
56 #include "AliRun.h"
57 #include "AliDigitizationInput.h"
58 #include "AliRunLoader.h"
59 #include "AliCDBManager.h"
60 #include "AliCDBEntry.h"
61 #include "AliEMCALDigit.h"
62 #include "AliEMCAL.h"
63 #include "AliEMCALLoader.h"
64 #include "AliEMCALDigitizer.h"
65 #include "AliEMCALSDigitizer.h"
66 #include "AliEMCALGeometry.h"
67 #include "AliEMCALCalibData.h"
68 #include "AliEMCALSimParam.h"
69 #include "AliEMCALRawDigit.h"
70 #include "AliCaloCalibPedestal.h"
71
72   namespace
73   {
74     Double_t HeavisideTheta(Double_t x)
75     {
76       Double_t signal = 0.;
77       
78       if (x > 0.) signal = 1.;  
79       
80       return signal;  
81     }
82     
83     Double_t AnalogFastORFunction(Double_t *x, Double_t *par)
84     {
85       Double_t v0 = par[0];
86       Double_t t0 = par[1];
87       Double_t tr = par[2];
88       
89       Double_t R1 = 1000.;
90       Double_t C1 = 33e-12;
91       Double_t R2 = 1800;
92       Double_t C2 = 22e-12;
93       
94       Double_t t  =   x[0];
95       
96       return (((0.8*(-((TMath::Power(C1,2)*C2*TMath::Power(TMath::E(),(-t + t0)/(C1*R1))*
97                         TMath::Power(R1,2)*R2)/(C1*R1 - C2*R2)) + 
98                      C1*C2*R1*R2*(1 - (C2*TMath::Power(TMath::E(),(-t + t0)/(C2*R2))*R2)/(-(C1*R1) + C2*R2)))*v0*
99                 HeavisideTheta(t - t0))/tr 
100                - (0.8*(C1*C2*R1*R2 - 
101                        (TMath::Power(C1,2)*C2*TMath::Power(TMath::E(),(-1.*t + t0 + 1.25*tr)/(C1*R1))*
102                         TMath::Power(R1,2)*R2)/(C1*R1 - C2*R2) + 
103                        (C1*TMath::Power(C2,2)*TMath::Power(TMath::E(),(-1.*t + t0 + 1.25*tr)/(C2*R2))*
104                         R1*TMath::Power(R2,2))/(C1*R1 - C2*R2))*v0*
105                   HeavisideTheta(t - t0 - 1.25*tr))/tr)/(C2*R1));
106     }
107   }
108
109 ClassImp(AliEMCALDigitizer)
110   
111   
112 //____________________________________________________________________________ 
113 AliEMCALDigitizer::AliEMCALDigitizer()
114   : AliDigitizer("",""),
115     fDefaultInit(kTRUE),
116     fDigitsInRun(0),
117     fInit(0),
118     fInput(0),
119     fInputFileNames(0x0),
120     fEventNames(0x0),
121     fDigitThreshold(0),
122     fMeanPhotonElectron(0),
123     fGainFluctuations(0),
124     fPinNoise(0),
125     fTimeNoise(0),
126     fTimeDelay(0),
127     fTimeResolutionPar0(0),
128     fTimeResolutionPar1(0),
129     fADCchannelEC(0),
130     fADCpedestalEC(0),
131     fADCchannelECDecal(0),
132     fTimeChannel(0),
133     fTimeChannelDecal(0),
134     fNADCEC(0),
135     fEventFolderName(""),
136     fFirstEvent(0),
137     fLastEvent(0),
138     fCalibData(0x0),
139     fSDigitizer(0x0)
140 {
141   // ctor
142   InitParameters() ; 
143   fDigInput = 0 ;                     // We work in the standalone mode
144 }
145
146 //____________________________________________________________________________ 
147 AliEMCALDigitizer::AliEMCALDigitizer(TString alirunFileName, TString eventFolderName)
148   : AliDigitizer("EMCALDigitizer", alirunFileName),
149     fDefaultInit(kFALSE),
150     fDigitsInRun(0),
151     fInit(0),
152     fInput(0),
153     fInputFileNames(0), 
154     fEventNames(0), 
155     fDigitThreshold(0),
156     fMeanPhotonElectron(0),
157     fGainFluctuations(0),
158     fPinNoise(0),
159     fTimeNoise(0),
160     fTimeDelay(0),
161     fTimeResolutionPar0(0),
162     fTimeResolutionPar1(0),
163     fADCchannelEC(0),
164     fADCpedestalEC(0),
165     fADCchannelECDecal(0),
166     fTimeChannel(0),
167     fTimeChannelDecal(0),
168     fNADCEC(0),
169     fEventFolderName(eventFolderName),
170     fFirstEvent(0),
171     fLastEvent(0),
172     fCalibData(0x0),
173     fSDigitizer(0x0)
174 {
175   // ctor
176   InitParameters() ; 
177   Init() ;
178   fDigInput = 0 ;                     // We work in the standalone mode
179 }
180
181 //____________________________________________________________________________ 
182 AliEMCALDigitizer::AliEMCALDigitizer(const AliEMCALDigitizer & d) 
183   : AliDigitizer(d.GetName(),d.GetTitle()),
184     fDefaultInit(d.fDefaultInit),
185     fDigitsInRun(d.fDigitsInRun),
186     fInit(d.fInit),
187     fInput(d.fInput),
188     fInputFileNames(d.fInputFileNames),
189     fEventNames(d.fEventNames),
190     fDigitThreshold(d.fDigitThreshold),
191     fMeanPhotonElectron(d.fMeanPhotonElectron),
192     fGainFluctuations(d.fGainFluctuations),
193     fPinNoise(d.fPinNoise),
194     fTimeNoise(d.fTimeNoise),
195     fTimeDelay(d.fTimeDelay),
196     fTimeResolutionPar0(d.fTimeResolutionPar0),
197     fTimeResolutionPar1(d.fTimeResolutionPar1),
198     fADCchannelEC(d.fADCchannelEC),
199     fADCpedestalEC(d.fADCpedestalEC),
200     fADCchannelECDecal(d.fADCchannelECDecal),
201     fTimeChannel(d.fTimeChannel), fTimeChannelDecal(d.fTimeChannelDecal),
202     fNADCEC(d.fNADCEC),
203     fEventFolderName(d.fEventFolderName),
204     fFirstEvent(d.fFirstEvent),
205     fLastEvent(d.fLastEvent),
206     fCalibData(d.fCalibData),
207     fSDigitizer(d.fSDigitizer ? new AliEMCALSDigitizer(*d.fSDigitizer) : 0)
208 {
209   // copyy ctor 
210 }
211
212 //____________________________________________________________________________ 
213 AliEMCALDigitizer::AliEMCALDigitizer(AliDigitizationInput * rd)
214   : AliDigitizer(rd,"EMCALDigitizer"),
215     fDefaultInit(kFALSE),
216     fDigitsInRun(0),
217     fInit(0),
218     fInput(0),
219     fInputFileNames(0),
220     fEventNames(0),
221     fDigitThreshold(0),
222     fMeanPhotonElectron(0),
223     fGainFluctuations(0),
224     fPinNoise(0.),
225     fTimeNoise(0.),
226     fTimeDelay(0.),
227     fTimeResolutionPar0(0.),
228     fTimeResolutionPar1(0.),
229     fADCchannelEC(0),
230     fADCpedestalEC(0),
231     fADCchannelECDecal(0),
232     fTimeChannel(0),
233     fTimeChannelDecal(0),    
234     fNADCEC(0),
235     fEventFolderName(0),
236     fFirstEvent(0),
237     fLastEvent(0),
238     fCalibData(0x0),
239     fSDigitizer(0x0)
240 {
241   // ctor Init() is called by RunDigitizer
242   fDigInput = rd ; 
243   fEventFolderName = fDigInput->GetInputFolderName(0) ;
244   SetTitle(dynamic_cast<AliStream*>(fDigInput->GetInputStream(0))->GetFileName(0));
245   InitParameters() ; 
246 }
247
248 //____________________________________________________________________________ 
249   AliEMCALDigitizer::~AliEMCALDigitizer()
250 {
251   //dtor
252   delete [] fInputFileNames ; 
253   delete [] fEventNames ; 
254   if (fSDigitizer) delete fSDigitizer;
255 }
256
257 //____________________________________________________________________________
258 void AliEMCALDigitizer::Digitize(Int_t event) 
259
260   
261   // Makes the digitization of the collected summable digits
262   // for this it first creates the array of all EMCAL modules
263   // filled with noise and after that adds contributions from 
264   // SDigits. This design helps to avoid scanning over the 
265   // list of digits to add  contribution of any new SDigit.
266   //
267   // JLK 26-Jun-2008
268   // Note that SDigit energy info is stored as an amplitude, so we
269   // must call the Calibrate() method of the SDigitizer to convert it
270   // back to an energy in GeV before adding it to the Digit
271   //
272   static int nEMC=0; //max number of digits possible
273   AliRunLoader *rl = AliRunLoader::Instance();
274   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
275   
276   if(emcalLoader){
277     Int_t readEvent = event ; 
278     // fDigInput is data member from AliDigitizer
279     if (fDigInput) 
280       readEvent = dynamic_cast<AliStream*>(fDigInput->GetInputStream(0))->GetCurrentEventNumber() ; 
281     AliDebug(1,Form("Adding event %d from input stream 0 %s %s", 
282                     readEvent, GetTitle(), fEventFolderName.Data())) ; 
283     rl->GetEvent(readEvent);
284     
285     TClonesArray * digits = emcalLoader->Digits() ; 
286     digits->Delete() ;  //correct way to clear array when memory is
287     //allocated by objects stored in array
288     
289     // Load Geometry
290     AliEMCALGeometry *geom = 0;
291     if (!rl->GetAliRun()) {
292       AliFatal("Could not get AliRun from runLoader");
293     }
294     else{
295       AliEMCAL * emcal  = (AliEMCAL*)rl->GetAliRun()->GetDetector("EMCAL");
296       geom = emcal->GetGeometry();
297       nEMC = geom->GetNCells();
298       AliDebug(1,Form("nEMC %i (number cells in EMCAL) | %s \n", nEMC, geom->GetName()));
299     }
300     
301     Int_t absID = -1 ;
302     
303     digits->Expand(nEMC) ;
304     
305     // RS create a digitizer on the fly
306     if (!fSDigitizer) fSDigitizer = new AliEMCALSDigitizer(rl->GetFileName().Data());
307     fSDigitizer->SetEventRange(0, -1) ;
308     //
309     //take all the inputs to add together and load the SDigits
310     TObjArray * sdigArray = new TObjArray(fInput) ;
311     sdigArray->AddAt(emcalLoader->SDigits(), 0) ;
312     Int_t i ;
313     Int_t embed = kFALSE;
314     for(i = 1 ; i < fInput ; i++){
315       TString tempo(fEventNames[i]) ; 
316       tempo += i ;
317         
318       AliRunLoader *  rl2 = AliRunLoader::GetRunLoader(tempo) ; 
319         
320       if (rl2==0) 
321         rl2 = AliRunLoader::Open(fInputFileNames[i], tempo) ; 
322         
323       if (fDigInput) 
324         readEvent = dynamic_cast<AliStream*>(fDigInput->GetInputStream(i))->GetCurrentEventNumber() ; 
325       Info("Digitize", "Adding event %d from input stream %d %s %s", readEvent, i, fInputFileNames[i].Data(), tempo.Data()) ; 
326       rl2->LoadSDigits();
327       //     rl2->LoadDigits();
328       rl2->GetEvent(readEvent);
329       if(rl2->GetDetectorLoader("EMCAL"))
330         {
331           AliEMCALLoader *emcalLoader2 = dynamic_cast<AliEMCALLoader*>(rl2->GetDetectorLoader("EMCAL"));
332  
333           if(emcalLoader2) {
334             //Merge 2 simulated sdigits
335             if(emcalLoader2->SDigits()){
336               TClonesArray* sdigits2 = emcalLoader2->SDigits();
337               sdigArray->AddAt(sdigits2, i) ;
338               // Check if first sdigit is of embedded type, if so, handle the sdigits differently:
339               // do not smear energy of embedded, do not add noise to any sdigits
340               if(sdigits2->GetEntriesFast()>0){
341                 //printf("Merged digit type: %d\n",dynamic_cast<AliEMCALDigit*> (sdigits2->At(0))->GetType());
342                 AliEMCALDigit * digit2 = dynamic_cast<AliEMCALDigit*> (sdigits2->At(0));
343                 if(digit2 && digit2->GetType()==AliEMCALDigit::kEmbedded){
344                   embed = kTRUE;
345                 }
346               }
347             }
348             
349           }//loader2
350           else  AliFatal("EMCAL Loader of second event not available!");
351           
352         }
353       else Fatal("Digitize", "Loader of second input not found");
354     }// input loop
355       
356     //Find the first tower with signal
357     Int_t nextSig = nEMC + 1 ; 
358     TClonesArray * sdigits ;  
359     for(i = 0 ; i < fInput ; i++){
360       if(i > 0 && embed) continue;
361       sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
362       if(sdigits){
363         if (sdigits->GetEntriesFast() ){
364           AliEMCALDigit *sd = dynamic_cast<AliEMCALDigit *>(sdigits->At(0));
365           if(sd){
366             Int_t curNext = sd->GetId() ;
367             if(curNext < nextSig) 
368               nextSig = curNext ;
369             AliDebug(1,Form("input %i : #sdigits %i \n",
370                             i, sdigits->GetEntriesFast()));
371               
372           }//First entry is not NULL
373           else {
374             Info("Digitize", "NULL sdigit pointer");
375             continue;
376           }
377         }//There is at least one entry
378         else {
379           Info("Digitize", "NULL sdigits array");
380           continue;
381         }
382       }// SDigits array exists
383       else {
384         Info("Digitizer","Null sdigit pointer");
385         continue;
386       }
387     }// input loop
388     AliDebug(1,Form("FIRST tower with signal %i \n", nextSig));
389       
390     TArrayI index(fInput) ;
391     index.Reset() ;  //Set all indexes to zero
392       
393     AliEMCALDigit * digit ;
394     AliEMCALDigit * curSDigit ;
395       
396     //Put Noise contribution, smear time and energy     
397     Float_t timeResolution = 0;
398     for(absID = 0; absID < nEMC; absID++){ // Nov 30, 2006 by PAI; was from 1 to nEMC
399         
400       if (IsDead(absID)) continue; // Don't instantiate dead digits
401         
402       Float_t energy = 0 ;    
403         
404       // amplitude set to zero, noise will be added later
405       Float_t noiseTime = 0.;
406       if(!embed) noiseTime = TimeOfNoise(); //No need for embedded events?
407       new((*digits)[absID]) AliEMCALDigit( -1, -1, absID, 0., noiseTime,kFALSE); // absID-1->absID
408       //look if we have to add signal?
409       digit = dynamic_cast<AliEMCALDigit *>(digits->At(absID)); // absID-1->absID
410         
411       if (digit) {
412           
413         if(absID==nextSig){
414
415           // Calculate time as time of the largest digit
416           Float_t time = digit->GetTime() ;
417           Float_t aTime= digit->GetAmplitude() ;
418             
419           // loop over input
420           Int_t nInputs = fInput;
421           if(embed) nInputs = 1; // In case of embedding, merge later real digits, do not smear energy and time of data
422           for(i = 0; i< nInputs ; i++){  //loop over (possible) merge sources
423             TClonesArray* sdtclarr = dynamic_cast<TClonesArray *>(sdigArray->At(i));
424             if(sdtclarr) {
425               Int_t sDigitEntries = sdtclarr->GetEntriesFast();
426                 
427               if(sDigitEntries > index[i] )
428                 curSDigit = dynamic_cast<AliEMCALDigit*>(sdtclarr->At(index[i])) ;      
429               else
430                 curSDigit = 0 ;
431                 
432               //May be several digits will contribute from the same input
433               while(curSDigit && (curSDigit->GetId() == absID)){           
434                 //Shift primary to separate primaries belonging different inputs
435                 Int_t primaryoffset ;
436                 if(fDigInput)
437                   primaryoffset = fDigInput->GetMask(i) ; 
438                 else
439                   primaryoffset = i ;
440                 curSDigit->ShiftPrimary(primaryoffset) ;
441                   
442                 if(curSDigit->GetAmplitude()>aTime) {
443                   aTime = curSDigit->GetAmplitude();
444                   time  = curSDigit->GetTime();
445                 }
446                   
447                 *digit = *digit + *curSDigit ;  //adds amplitudes of each digit
448                   
449                 index[i]++ ;
450                   
451                 if( sDigitEntries > index[i] )
452                   curSDigit = dynamic_cast<AliEMCALDigit*>(sdtclarr->At(index[i])) ;
453                 else
454                   curSDigit = 0 ;
455               }// while
456             }// source exists
457           }// loop over merging sources
458             
459             
460           //Here we convert the summed amplitude to an energy in GeV only for simulation or mixing of simulations
461           energy = fSDigitizer->Calibrate(digit->GetAmplitude()) ; // GeV
462             
463           // add fluctuations for photo-electron creation
464     // corrected fluctuations after comparison with beam test, Paraskevi Ganoti (06/11/2011)
465     Float_t fluct = static_cast<Float_t>((energy*fMeanPhotonElectron)/fGainFluctuations);
466     energy       *= static_cast<Float_t>(gRandom->Poisson(fluct)) / fluct ;    
467     
468           //calculate and set time
469           digit->SetTime(time) ;
470             
471           //Find next signal module
472           nextSig = nEMC + 1 ;
473           for(i = 0 ; i < nInputs ; i++){
474             sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
475               
476             if(sdigits){
477               Int_t curNext = nextSig ;
478               if(sdigits->GetEntriesFast() > index[i])
479                 {
480                   AliEMCALDigit * tmpdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(index[i]));
481                   if(tmpdigit)
482                     {
483                       curNext = tmpdigit->GetId() ;       
484                     }
485                 }
486               if(curNext < nextSig) nextSig = curNext ;
487             }// sdigits exist
488           } // input loop       
489             
490         }//absID==nextSig
491           
492         // add the noise now, no need for embedded events with real data
493         if(!embed)
494           energy += gRandom->Gaus(0., fPinNoise) ;
495           
496
497         // JLK 26-June-2008
498         //Now digitize the energy according to the fSDigitizer method,
499         //which merely converts the energy to an integer.  Later we will
500         //check that the stored value matches our allowed dynamic ranges
501         digit->SetAmplitude(fSDigitizer->Digitize(energy)) ; 
502           
503         //Set time resolution with final energy
504         timeResolution = GetTimeResolution(energy);
505         digit->SetTime(gRandom->Gaus(digit->GetTime(),timeResolution) ) ;
506         AliDebug(10,Form(" absID %5i energy %f nextSig %5i\n",
507                          absID, energy, nextSig));
508         //Add delay to time
509         digit->SetTime(digit->GetTime()+fTimeDelay) ;
510           
511       }// Digit pointer exists
512       else{
513         Info("Digitizer","Digit pointer is null");
514       }
515     } // for(absID = 0; absID < nEMC; absID++)
516       
517       //ticks->Delete() ;
518       //delete ticks ;
519
520       //Embed simulated amplitude (and time?) to real data digits
521     if(embed){
522       for(Int_t input = 1; input<fInput; input++){
523         TClonesArray *realDigits = dynamic_cast<TClonesArray*> (sdigArray->At(input));
524         if(!realDigits){
525           AliDebug(1,"No sdigits to merge\n");
526           continue;
527         }
528         for(Int_t i2 = 0 ; i2 < realDigits->GetEntriesFast() ; i2++){
529           AliEMCALDigit * digit2 = dynamic_cast<AliEMCALDigit*>( realDigits->At(i2) ) ;
530           if(digit2){
531             digit = dynamic_cast<AliEMCALDigit*>( digits->At(digit2->GetId()) ) ;
532             if(digit){
533                 
534               // Put the embedded cell energy in same units as simulated sdigits -> transform to energy units
535               // and multiply to get a big int.
536               Float_t time2 = digit2->GetTime(); 
537               Float_t e2    = digit2->GetAmplitude();
538               CalibrateADCTime(e2,time2,digit2->GetId());
539               Float_t amp2 = fSDigitizer->Digitize(e2);
540                 
541               Float_t e0    = digit ->GetAmplitude();
542               if(e0>0.01)
543                 AliDebug(1,Form("digit 1: Abs ID %d, amp %f, type %d, time %e; digit2: Abs ID %d, amp %f, type %d, time %e\n",
544                                 digit ->GetId(),digit ->GetAmplitude(), digit ->GetType(), digit->GetTime(),
545                                 digit2->GetId(),amp2,                   digit2->GetType(), time2           ));
546                 
547               // Sum amplitudes, change digit amplitude
548               digit->SetAmplitude( digit->GetAmplitude() + amp2);
549               // Rough assumption, assign time of the largest of the 2 digits, 
550               // data or signal to the final digit.
551               if(amp2 > digit->GetAmplitude())  digit->SetTime(time2);
552               //Mark the new digit as embedded
553               digit->SetType(AliEMCALDigit::kEmbedded);
554                 
555               if(digit2->GetAmplitude()>0.01 && e0> 0.01 )
556                 AliDebug(1,Form("Embedded digit: Abs ID %d, amp %f, type %d\n",
557                                 digit->GetId(), digit->GetAmplitude(), digit->GetType()));
558             }//digit2 
559           }//digit2
560         }//loop digit 2
561       }//input loop
562     }//embed
563       
564     //JLK is it better to call Clear() here?
565     delete sdigArray ; //We should not delete its contents
566       
567     //remove digits below thresholds
568     // until 10-02-2010 remove digits with energy smaller than fDigitThreshold 3*fPinNoise
569     // now, remove digits with Digitized ADC smaller than fDigitThreshold = 3,
570     // before merge in the same loop real data digits if available
571     Float_t energy = 0;
572     Float_t time   = 0;
573     for(i = 0 ; i < nEMC ; i++){
574       digit = dynamic_cast<AliEMCALDigit*>( digits->At(i) ) ;
575       if(digit){
576           
577         //Then get the energy in GeV units.
578         energy = fSDigitizer->Calibrate(digit->GetAmplitude()) ;
579         //Then digitize using the calibration constants of the ocdb
580         Float_t ampADC = energy;
581         DigitizeEnergyTime(ampADC, time, digit->GetId())  ;       
582         if(ampADC < fDigitThreshold)
583           digits->RemoveAt(i) ;
584           
585       }// digit exists
586     } // digit loop
587       
588     digits->Compress() ;  
589       
590     Int_t ndigits = digits->GetEntriesFast() ; 
591
592     //JLK 26-June-2008
593     //After we have done the summing and digitizing to create the
594     //digits, now we want to calibrate the resulting amplitude to match
595     //the dynamic range of our real data.  
596     for (i = 0 ; i < ndigits ; i++) { 
597       digit = dynamic_cast<AliEMCALDigit *>( digits->At(i) ) ; 
598       if(digit){
599         digit->SetIndexInList(i) ; 
600         energy = fSDigitizer->Calibrate(digit->GetAmplitude()) ;
601         time   = digit->GetTime();
602         Float_t ampADC = energy;
603         DigitizeEnergyTime(ampADC, time, digit->GetId());
604         digit->SetAmplitude(ampADC) ;
605         digit->SetTime(time);
606         // printf("digit amplitude set at end: i %d, amp %f\n",i,digit->GetAmplitude());
607       }// digit exists
608     }//Digit loop
609       
610   }
611   else AliFatal("EMCALLoader is NULL!") ;
612   
613 }
614
615 //_____________________________________________________________________
616 void AliEMCALDigitizer::DigitizeEnergyTime(Float_t & energy, Float_t & time, const Int_t absId)
617
618   // JLK 26-June-2008
619   // Returns digitized value of the energy in a cell absId
620   // using the calibration constants stored in the OCDB
621   // or default values if no CalibData object is found.
622   // This effectively converts everything to match the dynamic range
623   // of the real data we will collect
624   //
625   // Load Geometry
626   const AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
627   
628   if (geom==0){
629     AliFatal("Did not get geometry from EMCALLoader");
630     return;
631   }
632   
633   Int_t iSupMod = -1;
634   Int_t nModule = -1;
635   Int_t nIphi   = -1;
636   Int_t nIeta   = -1;
637   Int_t iphi    = -1;
638   Int_t ieta    = -1;
639   Bool_t bCell = geom->GetCellIndex(absId, iSupMod, nModule, nIphi, nIeta) ;
640   if(!bCell)
641     Error("DigitizeEnergyTime","Wrong cell id number : absId %i ", absId) ;
642   geom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
643   
644   if(fCalibData) {
645     fADCpedestalEC     = fCalibData->GetADCpedestal     (iSupMod,ieta,iphi);
646     fADCchannelEC      = fCalibData->GetADCchannel      (iSupMod,ieta,iphi);
647     fADCchannelECDecal = fCalibData->GetADCchannelDecal (iSupMod,ieta,iphi);
648     fTimeChannel       = fCalibData->GetTimeChannel     (iSupMod,ieta,iphi,0); // Assign bunch crossing number equal to 0 (has simulation different bunch crossings?)
649     fTimeChannelDecal  = fCalibData->GetTimeChannelDecal(iSupMod,ieta,iphi);      
650   }
651   
652   //Apply calibration to get ADC counts and partial decalibration as especified in OCDB
653   energy = (energy + fADCpedestalEC)/fADCchannelEC/fADCchannelECDecal   ;
654   time  += fTimeChannel-fTimeChannelDecal;
655   
656   if(energy > fNADCEC ) 
657     energy =  fNADCEC ; 
658 }
659
660 //_____________________________________________________________________
661 void AliEMCALDigitizer::Decalibrate(AliEMCALDigit *digit)
662
663         // Load Geometry
664         const AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
665         
666         if (!geom) {
667                 AliFatal("Did not get geometry from EMCALLoader");
668                 return;
669         }
670         
671         Int_t absId   = digit->GetId();
672         Int_t iSupMod = -1;
673         Int_t nModule = -1;
674         Int_t nIphi   = -1;
675         Int_t nIeta   = -1;
676         Int_t iphi    = -1;
677         Int_t ieta    = -1;
678         
679         Bool_t bCell = geom->GetCellIndex(absId, iSupMod, nModule, nIphi, nIeta) ;
680         
681         if (!bCell) Error("Decalibrate","Wrong cell id number : absId %i ", absId) ;
682         geom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
683         
684         if (fCalibData) {
685                 fADCchannelEC = fCalibData->GetADCchannel(iSupMod,ieta,iphi);
686                 float factor = 1./(fADCchannelEC/0.0162);
687                 
688                 *digit = *digit * factor;
689         }
690 }
691
692 //_____________________________________________________________________
693 void AliEMCALDigitizer::CalibrateADCTime(Float_t & adc, Float_t & time, const Int_t absId)
694
695   // Returns the energy in a cell absId with a given adc value
696   // using the calibration constants stored in the OCDB. Time also corrected from parameter in OCDB
697   // Used in case of embedding, transform ADC counts from real event into energy
698   // so that we can add the energy of the simulated sdigits which are in energy
699   // units.
700   // Same as in AliEMCALClusterizer::Calibrate() but here we do not reject channels being marked as hot
701   // or with time out of window
702   
703   // Load Geometry
704   const AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
705   
706   if (!geom){
707     AliFatal("Did not get geometry from EMCALLoader");
708     return;
709   }
710   
711   Int_t iSupMod = -1;
712   Int_t nModule = -1;
713   Int_t nIphi   = -1;
714   Int_t nIeta   = -1;
715   Int_t iphi    = -1;
716   Int_t ieta    = -1;
717   Bool_t bCell = geom->GetCellIndex(absId, iSupMod, nModule, nIphi, nIeta) ;
718   if(!bCell)
719     Error("CalibrateADCTime","Wrong cell id number : absId %i ", absId) ;
720   geom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
721   
722   if(fCalibData) {
723     fADCpedestalEC = fCalibData->GetADCpedestal(iSupMod,ieta,iphi);
724     fADCchannelEC  = fCalibData->GetADCchannel (iSupMod,ieta,iphi);
725     fTimeChannel   = fCalibData->GetTimeChannel(iSupMod,ieta,iphi,0);// Assign bunch crossing number equal to 0 (has simulation different bunch crossings?)
726   }
727   
728   adc   = adc * fADCchannelEC - fADCpedestalEC;
729   time -= fTimeChannel;
730   
731   
732 }
733
734
735 //____________________________________________________________________________
736 void AliEMCALDigitizer::Digitize(Option_t *option) 
737
738   // Steering method to process digitization for events
739   // in the range from fFirstEvent to fLastEvent.
740   // This range is optionally set by SetEventRange().
741   // if fLastEvent=-1, then process events until the end.
742   // by default fLastEvent = fFirstEvent (process only one event)
743
744   if (!fInit) { // to prevent overwrite existing file
745     Error( "Digitize", "Give a version name different from %s", fEventFolderName.Data() ) ;
746     return ;
747   } 
748
749   if (strstr(option,"print")) {
750
751     Print();
752     return ; 
753   }
754   
755   if(strstr(option,"tim"))
756     gBenchmark->Start("EMCALDigitizer");
757
758   AliRunLoader *rl = AliRunLoader::Instance();
759   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
760   Int_t nEvents = 0;
761   if(!emcalLoader){
762     AliFatal("Did not get the  Loader");
763   }
764   else{
765     
766     if (fLastEvent == -1)  {
767       fLastEvent = rl->GetNumberOfEvents() - 1 ;
768     }
769     else if (fDigInput) 
770       fLastEvent = fFirstEvent ; // what is this ??
771     
772     nEvents = fLastEvent - fFirstEvent + 1;
773     Int_t ievent  = -1;
774
775     TClonesArray* digitsTMP = new TClonesArray("AliEMCALDigit",    32*96);
776     TClonesArray* digitsTRG = new TClonesArray("AliEMCALRawDigit", 32*96);
777
778     rl->LoadSDigits("EMCAL");
779     for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
780       
781       rl->GetEvent(ievent);
782       
783       Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
784       
785       WriteDigits() ;
786       
787       //Trigger Digits
788       //-------------------------------------
789      
790
791       Digits2FastOR(digitsTMP, digitsTRG);  
792       
793       WriteDigits(digitsTRG);
794       
795       (emcalLoader->TreeD())->Fill();
796       
797       emcalLoader->WriteDigits(   "OVERWRITE");
798       
799       Unload();
800       
801       digitsTRG  ->Delete();
802       digitsTMP  ->Delete();
803       
804       //-------------------------------------
805       
806       if(strstr(option,"deb"))
807         PrintDigits(option);
808       if(strstr(option,"table")) gObjectTable->Print();
809       
810       //increment the total number of Digits per run 
811       fDigitsInRun += emcalLoader->Digits()->GetEntriesFast() ;  
812     }//loop
813        
814   }//loader exists
815   
816   if(strstr(option,"tim")){
817     gBenchmark->Stop("EMCALDigitizer");
818     Float_t cputime   = gBenchmark->GetCpuTime("EMCALDigitizer");
819     Float_t avcputime = cputime;
820     if(nEvents==0) avcputime = 0 ;
821     AliInfo(Form("Digitize: took %f seconds for Digitizing %f seconds per event", cputime, avcputime)) ;
822   } 
823 }
824
825 //__________________________________________________________________
826 Float_t AliEMCALDigitizer::GetTimeResolution(Float_t energy) const
827 {  
828   // From F. Blanco
829   Float_t res = -1;
830   if (energy > 0) {
831     res = TMath::Sqrt(fTimeResolutionPar0 + 
832                       fTimeResolutionPar1/(energy*energy) );
833   }
834   // parametrization above is for ns. Convert to seconds:
835   return res*1e-9;
836 }
837
838 //____________________________________________________________________________ 
839 void AliEMCALDigitizer::Digits2FastOR(TClonesArray* digitsTMP, TClonesArray* digitsTRG)
840 {
841   // FEE digits afterburner to produce TRG digits 
842   // we are only interested in the FEE digit deposited energy
843   // to be converted later into a voltage value
844   
845   // push the FEE digit to its associated FastOR (numbered from 0:95)
846   // TRU is in charge of summing module digits
847   
848   AliRunLoader *runLoader = AliRunLoader::Instance();
849   
850   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(runLoader->GetDetectorLoader("EMCAL"));
851   if (!emcalLoader) AliFatal("Did not get the  Loader");
852   
853     const  AliEMCALGeometry* geom = AliEMCALGeometry::GetInstance();
854       
855       // build FOR from simulated digits
856       // and xfer to the corresponding TRU input (mapping)
857       
858           TClonesArray* sdigits = emcalLoader->SDigits();
859                 
860           TClonesArray *digits = (TClonesArray*)sdigits->Clone();
861                 
862           AliDebug(999,Form("=== %d SDigits to trigger digits ===",digits->GetEntriesFast()));
863                 
864       TIter NextDigit(digits);
865       while (AliEMCALDigit* digit = (AliEMCALDigit*)NextDigit())
866       {
867         if (IsDead(digit)) continue;
868         
869               Decalibrate(digit);   
870                   
871         Int_t id = digit->GetId();
872         
873         Int_t trgid;
874         if (geom && geom->GetFastORIndexFromCellIndex(id , trgid)) 
875         {
876           AliDebug(1,Form("trigger digit id: %d from cell id: %d\n",trgid,id));
877           
878           AliEMCALDigit* d = static_cast<AliEMCALDigit*>(digitsTMP->At(trgid));
879           
880           if (!d)
881           {
882             new((*digitsTMP)[trgid]) AliEMCALDigit(*digit);
883             d = (AliEMCALDigit*)digitsTMP->At(trgid);
884             d->SetId(trgid);
885           }     
886           else
887           {
888             *d = *d + *digit;
889           }
890         }
891       }
892       
893       if (AliDebugLevel()) printf("Number of TRG digits: %d\n",digitsTMP->GetEntriesFast());
894                 
895       Int_t    nSamples = 32;
896       Int_t *timeSamples = new Int_t[nSamples];
897       
898       NextDigit = TIter(digitsTMP);
899       while (AliEMCALDigit* digit = (AliEMCALDigit*)NextDigit())
900       {
901         if (digit)
902         {
903           Int_t     id = digit->GetId();
904           Float_t time = 50.e-9;
905           
906           Double_t depositedEnergy = 0.;
907           for (Int_t j = 1; j <= digit->GetNprimary(); j++) depositedEnergy += digit->GetDEPrimary(j);
908           
909           if (AliDebugLevel()) printf("Deposited Energy: %f\n", depositedEnergy);
910           
911           // FIXME: Check digit time!
912           if (depositedEnergy)
913           {
914                         depositedEnergy += gRandom->Gaus(0., .08);
915                           
916             DigitalFastOR(time, depositedEnergy, timeSamples, nSamples);
917             
918             for (Int_t j=0;j<nSamples;j++) 
919             {
920               timeSamples[j] = ((j << 12) & 0xFF000) | (timeSamples[j] & 0xFFF);
921             }
922             
923             new((*digitsTRG)[digitsTRG->GetEntriesFast()]) AliEMCALRawDigit(id, timeSamples, nSamples);
924                           
925                         if (AliDebugLevel()) ((AliEMCALRawDigit*)digitsTRG->At(digitsTRG->GetEntriesFast() - 1))->Print("");
926           }
927         }
928       }
929       
930       delete [] timeSamples;
931                 
932           if (digits && digits->GetEntriesFast()) digits->Delete();
933 }
934
935 //____________________________________________________________________________ 
936 void AliEMCALDigitizer::DigitalFastOR( Double_t time, Double_t dE, Int_t timeSamples[], Int_t nSamples )
937 {
938         // parameters:  
939         // id: 0..95
940         const Int_t    reso = 11;      // 11-bit resolution ADC
941         const Double_t vFSR = 1;       // Full scale input voltage range
942 //      const Double_t dNe  = 125;     // signal of the APD per MeV of energy deposit in a tower: 125 photo-e-/MeV @ M=30
943         const Double_t dNe  = 125/1.3; // F-ALTRO max V. FEE: factor 4 
944         const Double_t vA   = .136e-6; // CSP output range: 0.136uV/e-
945         const Double_t rise = 50e-9;   // rise time (10-90%) of the FastOR signal before shaping
946         
947         const Double_t kTimeBinWidth = 25E-9; // sampling frequency (40MHz)
948         
949         Double_t vV = 1000. * dE * dNe * vA; // GeV 2 MeV
950                 
951         TF1 signalF("signal", AnalogFastORFunction, 0, nSamples * kTimeBinWidth, 3);
952         signalF.SetParameter( 0,   vV ); 
953         signalF.SetParameter( 1, time ); // FIXME: when does the signal arrive? Might account for cable lengths
954         signalF.SetParameter( 2, rise );
955         
956         for (Int_t iTime=0; iTime<nSamples; iTime++) 
957         {
958                 // FIXME: add noise (probably not simply Gaussian) according to DA measurements
959                 // probably plan an access to OCDB
960                 
961                 timeSamples[iTime] = int((TMath::Power(2, reso) / vFSR) * signalF.Eval(iTime * kTimeBinWidth) + 0.5);
962         }
963 }
964
965 //____________________________________________________________________________ 
966 Bool_t AliEMCALDigitizer::Init()
967 {
968   // Makes all memory allocations
969   fInit = kTRUE ; 
970   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
971   
972   if ( emcalLoader == 0 ) {
973     Fatal("Init", "Could not obtain the AliEMCALLoader");  
974     return kFALSE;
975   } 
976   
977   fFirstEvent = 0 ; 
978   fLastEvent = fFirstEvent ; 
979   
980   if(fDigInput)
981     fInput = fDigInput->GetNinputs() ; 
982   else 
983     fInput           = 1 ;
984
985   fInputFileNames  = new TString[fInput] ; 
986   fEventNames      = new TString[fInput] ; 
987   fInputFileNames[0] = GetTitle() ; 
988   fEventNames[0]     = fEventFolderName.Data() ; 
989   Int_t index ; 
990   for (index = 1 ; index < fInput ; index++) {
991     fInputFileNames[index] = dynamic_cast<AliStream*>(fDigInput->GetInputStream(index))->GetFileName(0); 
992     TString tempo = fDigInput->GetInputFolderName(index) ;
993     fEventNames[index] = tempo.Remove(tempo.Length()-1) ; // strip of the stream number added bt fDigInput 
994   }
995     
996   //Calibration instance
997   fCalibData = emcalLoader->CalibData();
998   return fInit ;    
999 }
1000
1001 //____________________________________________________________________________ 
1002 void AliEMCALDigitizer::InitParameters()
1003
1004   // Parameter initialization for digitizer
1005   
1006   // Get the parameters from the OCDB via the loader
1007   AliRunLoader *rl = AliRunLoader::Instance();
1008   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
1009   AliEMCALSimParam * simParam = 0x0;
1010   if(emcalLoader) simParam = emcalLoader->SimulationParameters();
1011         
1012   if(!simParam){
1013     simParam = AliEMCALSimParam::GetInstance();
1014     AliWarning("Simulation Parameters not available in OCDB?");
1015   }
1016   
1017   fMeanPhotonElectron = simParam->GetMeanPhotonElectron() ; // 4400;  // electrons per GeV 
1018   fGainFluctuations   = simParam->GetGainFluctuations()   ; // 15.0; 
1019   
1020   fPinNoise           = simParam->GetPinNoise();//0.012; // pin noise in GeV from analysis test beam data 
1021   if (fPinNoise < 0.0001 ) 
1022     Warning("InitParameters", "No noise added\n") ; 
1023   fTimeNoise          = simParam->GetTimeNoise(); // 1.28E-5 s
1024   fDigitThreshold     = simParam->GetDigitThreshold(); //fPinNoise * 3; // 3 * sigma
1025   fTimeResolutionPar0 = simParam->GetTimeResolutionPar0(); 
1026   fTimeResolutionPar1 = simParam->GetTimeResolutionPar1(); 
1027   fTimeDelay          = simParam->GetTimeDelay(); //600e-9 ; // 600 ns
1028   
1029   // These defaults are normally not used. 
1030   // Values are read from calibration database instead
1031   fADCchannelEC       = 0.0153; // Update 24 Apr 2007: 250./16/1024 - width of one ADC channel in GeV
1032   fADCpedestalEC      = 0.0 ;   // GeV
1033   fADCchannelECDecal  = 1.0;    // No decalibration by default
1034   fTimeChannel        = 0.0;    // No time calibration by default
1035   fTimeChannelDecal   = 0.0;    // No time decalibration by default
1036
1037   fNADCEC             = simParam->GetNADCEC();//(Int_t) TMath::Power(2,16) ;  // number of channels in Tower ADC - 65536
1038   
1039   AliDebug(2,Form("Mean Photon Electron %d, Gain Fluct. %2.1f; Noise: APD %f, Time %f; Digit Threshold %d,Time Resolution Par0 %g Par1 %g,NADCEC %d",
1040                   fMeanPhotonElectron, fGainFluctuations, fPinNoise,fTimeNoise, fDigitThreshold,fTimeResolutionPar0,fTimeResolutionPar1,fNADCEC));
1041   
1042 }
1043
1044 //__________________________________________________________________
1045 void AliEMCALDigitizer::Print1(Option_t * option)
1046 { // 19-nov-04 - just for convenience
1047   Print(); 
1048   PrintDigits(option);
1049 }
1050
1051 //__________________________________________________________________
1052 void AliEMCALDigitizer::Print (Option_t * ) const
1053 {
1054   // Print Digitizer's parameters
1055   printf("Print: \n------------------- %s -------------", GetName() ) ; 
1056   if( strcmp(fEventFolderName.Data(), "") != 0 ){
1057     printf(" Writing Digits to branch with title  %s\n", fEventFolderName.Data()) ;
1058     
1059     Int_t nStreams ; 
1060     if (fDigInput) 
1061       nStreams =  GetNInputStreams() ;
1062     else 
1063       nStreams = fInput ; 
1064     
1065     AliRunLoader *rl=0;
1066     
1067     Int_t index = 0 ;  
1068     for (index = 0 ; index < nStreams ; index++) {  
1069       TString tempo(fEventNames[index]) ; 
1070       tempo += index ;
1071       if ((rl = AliRunLoader::GetRunLoader(tempo)) == 0)
1072         rl = AliRunLoader::Open(fInputFileNames[index], tempo) ; 
1073       AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
1074       if(emcalLoader){
1075         TString fileName( emcalLoader->GetSDigitsFileName() ) ; 
1076         if ( fEventNames[index] != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name 
1077           fileName = fileName.ReplaceAll(".root", "") + "_" + fEventNames[index]  + ".root" ;
1078         printf ("Adding SDigits from %s %s\n", fInputFileNames[index].Data(), fileName.Data()) ; 
1079       }//loader
1080     }
1081     
1082     AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1083     
1084     if(emcalLoader) printf("\nWriting digits to %s", emcalLoader->GetDigitsFileName().Data()) ;
1085     else printf("\nNULL LOADER");
1086     
1087     printf("\nWith following parameters:\n") ;
1088     printf("    Electronics noise in EMC, APD (fPinNoise) = %f, Time = %f \n", fPinNoise, fTimeNoise) ;
1089     printf("    Threshold  in Tower  (fDigitThreshold) = %d\n", fDigitThreshold)  ;
1090     printf("---------------------------------------------------\n")  ;
1091   }
1092   else
1093     printf("Print: AliEMCALDigitizer not initialized") ; 
1094 }
1095
1096 //__________________________________________________________________
1097 void AliEMCALDigitizer::PrintDigits(Option_t * option)
1098 {
1099   //utility method for printing digit information
1100   
1101   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1102   if(emcalLoader){
1103     TClonesArray * digits  = emcalLoader->Digits() ;
1104     TClonesArray * sdigits = emcalLoader->SDigits() ;
1105     
1106     printf("\n #Digits: %d : sdigits %d ", digits->GetEntriesFast(), sdigits->GetEntriesFast()) ; 
1107     printf("\n event %d", emcalLoader->GetRunLoader()->GetEventNumber());
1108     
1109     if(strstr(option,"all")){  
1110       //loop over digits
1111       AliEMCALDigit * digit;
1112       printf("\nEMC digits (with primaries):\n")  ;
1113       printf("\n   Id  Amplitude    Time          Index Nprim: Primaries list \n") ;    
1114       Int_t index ;
1115       for (index = 0 ; index < digits->GetEntries() ; index++) {
1116         digit = dynamic_cast<AliEMCALDigit *>(digits->At(index)) ;
1117         if(digit){
1118           printf("\n%6d  %8f    %6.5e %4d      %2d : ",
1119                  digit->GetId(), digit->GetAmplitude(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;  
1120           Int_t iprimary;
1121           for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
1122             printf("%d ",digit->GetPrimary(iprimary+1) ) ; 
1123           }
1124         }// digit exists
1125       }// loop
1126     }
1127     printf("\n");
1128   }// loader exists
1129   else printf("NULL LOADER, cannot print\n");
1130 }
1131
1132 //__________________________________________________________________
1133 Float_t AliEMCALDigitizer::TimeOfNoise(void)
1134 {  
1135   // Calculates the time signal generated by noise
1136   //printf("Time noise %e\n",fTimeNoise);
1137   return gRandom->Rndm() * fTimeNoise;
1138 }
1139
1140 //__________________________________________________________________
1141 void AliEMCALDigitizer::Unload() 
1142 {  
1143   // Unloads the SDigits and Digits
1144   AliRunLoader *rl=0;
1145     
1146   Int_t i ; 
1147   for(i = 1 ; i < fInput ; i++){
1148     TString tempo(fEventNames[i]) ; 
1149     tempo += i;
1150     if ((rl = AliRunLoader::GetRunLoader(tempo))) 
1151       rl->GetDetectorLoader("EMCAL")->UnloadSDigits() ; 
1152   }
1153   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1154   if(emcalLoader)emcalLoader->UnloadDigits() ; 
1155   else AliFatal("Did not get the loader");
1156 }
1157
1158 //_________________________________________________________________________________________
1159 void AliEMCALDigitizer::WriteDigits()
1160 { // Makes TreeD in the output file. 
1161   // Check if branch already exists: 
1162   //   if yes, exit without writing: ROOT TTree does not support overwriting/updating of 
1163   //      already existing branches. 
1164   //   else creates branch with Digits, named "EMCAL", title "...",
1165   //      and branch "AliEMCALDigitizer", with the same title to keep all the parameters
1166   //      and names of files, from which digits are made.
1167   
1168   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1169   
1170   if(emcalLoader){
1171     const TClonesArray * digits = emcalLoader->Digits() ; 
1172     TTree * treeD = emcalLoader->TreeD(); 
1173     if ( !treeD ) {
1174       emcalLoader->MakeDigitsContainer();
1175       treeD = emcalLoader->TreeD(); 
1176     }
1177     
1178     // -- create Digits branch
1179     Int_t bufferSize = 32000 ;    
1180     TBranch * digitsBranch = 0;
1181     if ((digitsBranch = treeD->GetBranch("EMCAL"))) {
1182       digitsBranch->SetAddress(&digits);
1183       AliWarning("Digits Branch already exists. Not all digits will be visible");
1184     }
1185     else
1186       treeD->Branch("EMCAL","TClonesArray",&digits,bufferSize);
1187     //digitsBranch->SetTitle(fEventFolderName);
1188     
1189     //  treeD->Fill() ;
1190     /*  
1191      emcalLoader->WriteDigits("OVERWRITE");
1192      emcalLoader->WriteDigitizer("OVERWRITE");
1193      
1194      Unload() ; 
1195      */
1196     
1197   }// loader exists
1198   else AliFatal("Loader not available");
1199 }
1200
1201 //__________________________________________________________________
1202 void AliEMCALDigitizer::WriteDigits(TClonesArray* digits, const char* branchName)
1203 { // overloaded method
1204   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1205   if(emcalLoader){
1206     
1207     TTree* treeD = emcalLoader->TreeD(); 
1208     if (!treeD) 
1209       {
1210         emcalLoader->MakeDigitsContainer();
1211         treeD = emcalLoader->TreeD(); 
1212       }
1213     
1214     // -- create Digits branch
1215     Int_t bufferSize = 32000;
1216     
1217     if (TBranch* triggerBranch = treeD->GetBranch(branchName)) 
1218       {
1219         triggerBranch->SetAddress(&digits);
1220       }
1221     else
1222       {
1223         treeD->Branch(branchName,"TClonesArray",&digits,bufferSize);
1224       }
1225     
1226     //  treeD->Fill();
1227   }// loader exists
1228   else AliFatal("Loader not available");
1229 }
1230
1231 //__________________________________________________________________
1232 Bool_t AliEMCALDigitizer::IsDead(AliEMCALDigit *digit) 
1233 {
1234   AliRunLoader *runLoader = AliRunLoader::Instance();
1235   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(runLoader->GetDetectorLoader("EMCAL"));
1236   if (!emcalLoader) AliFatal("Did not get the  Loader");
1237   
1238   AliCaloCalibPedestal *caloPed = emcalLoader->PedestalData();
1239   if (!caloPed) {
1240     AliWarning("Could not access pedestal data! No dead channel removal applied");
1241     return kFALSE;
1242   }
1243   
1244         // Load Geometry
1245   const AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
1246   if (!geom) AliFatal("Did not get geometry from EMCALLoader");
1247   
1248   Int_t absId   = digit->GetId();
1249   Int_t iSupMod = -1;
1250   Int_t nModule = -1;
1251   Int_t nIphi   = -1;
1252   Int_t nIeta   = -1;
1253   Int_t iphi    = -1;
1254   Int_t ieta    = -1;
1255         
1256   Bool_t bCell = geom->GetCellIndex(absId, iSupMod, nModule, nIphi, nIeta) ;
1257         
1258   if (!bCell) Error("IsDead","Wrong cell id number : absId %i ", absId) ;
1259   geom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
1260
1261   Int_t channelStatus = (Int_t)(caloPed->GetDeadMap(iSupMod))->GetBinContent(ieta,iphi);
1262
1263   if (channelStatus == AliCaloCalibPedestal::kDead) 
1264     return kTRUE;
1265   else
1266     return kFALSE;
1267 }
1268
1269
1270 //__________________________________________________________________
1271 Bool_t AliEMCALDigitizer::IsDead(Int_t absId)
1272 {
1273     AliRunLoader *runLoader = AliRunLoader::Instance();
1274     AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(runLoader->GetDetectorLoader("EMCAL"));
1275     if (!emcalLoader) AliFatal("Did not get the  Loader");
1276     
1277     AliCaloCalibPedestal *caloPed = emcalLoader->PedestalData();
1278     if (!caloPed) {
1279         AliWarning("Could not access pedestal data! No dead channel removal applied");
1280         return kFALSE;
1281     }
1282     
1283         // Load Geometry
1284     const AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance();
1285     if (!geom) AliFatal("Did not get geometry from EMCALLoader");
1286     
1287     Int_t iSupMod = -1;
1288     Int_t nModule = -1;
1289     Int_t nIphi   = -1;
1290     Int_t nIeta   = -1;
1291     Int_t iphi    = -1;
1292     Int_t ieta    = -1;
1293         
1294     Bool_t bCell = geom->GetCellIndex(absId, iSupMod, nModule, nIphi, nIeta) ;
1295         
1296     if (!bCell) Error("IsDead","Wrong cell id number : absId %i ", absId) ;
1297     geom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
1298     
1299     Int_t channelStatus = (Int_t)(caloPed->GetDeadMap(iSupMod))->GetBinContent(ieta,iphi);
1300     
1301     if (channelStatus == AliCaloCalibPedestal::kDead)
1302         return kTRUE;
1303     else
1304         return kFALSE;
1305 }
1306