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