]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALDigitizer.cxx
Fixed static build, update rulechecker (Johny)
[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     AliFatal("Could not get AliRun from runLoader");
307   }
308   else{
309     AliEMCAL * emcal  = (AliEMCAL*)rl->GetAliRun()->GetDetector("EMCAL");
310     geom = emcal->GetGeometry();
311     nEMC = geom->GetNCells();
312     AliDebug(1,Form("nEMC %i (number cells in EMCAL) | %s \n", nEMC, geom->GetName()));
313   }
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         if(emcalLoader2) sdigArray->AddAt(emcalLoader2->SDigits(), i) ;
353         else             AliFatal("EMCAL Loader of second event not available!");
354       }
355       else Fatal("Digitize", "Loader of second input not found");
356     }
357     
358     //Find the first tower with signal
359     Int_t nextSig = nEMC + 1 ; 
360     TClonesArray * sdigits ;  
361     for(i = 0 ; i < fInput ; i++){
362       
363       sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
364       if(sdigits){
365         if (sdigits->GetEntriesFast() ){
366           if(dynamic_cast<AliEMCALDigit *>(sdigits->At(0))){
367             
368             Int_t curNext = dynamic_cast<AliEMCALDigit *>(sdigits->At(0))->GetId() ;
369             if(curNext < nextSig) 
370               nextSig = curNext ;
371             AliDebug(1,Form("input %i : #sdigits %i \n",
372                             i, sdigits->GetEntriesFast()));
373             
374           }//First entry is not NULL
375           else {
376             Info("Digitize", "NULL sdigit pointer");
377             continue;
378           }
379         }//There is at least one entry
380         else {
381           Info("Digitize", "NULL sdigits array");
382           continue;
383         }
384       }// SDigits array exists
385       else {
386         Info("Digitizer","Null sdigit pointer");
387         continue;
388       }
389     }// input loop
390     AliDebug(1,Form("FIRST tower with signal %i \n", nextSig));
391     
392     TArrayI index(fInput) ;
393     index.Reset() ;  //Set all indexes to zero
394     
395     AliEMCALDigit * digit ;
396     AliEMCALDigit * curSDigit ;
397     
398     //  TClonesArray * ticks = new TClonesArray("AliEMCALTick",1000) ;
399     
400     //Put Noise contribution
401     for(absID = 0; absID < nEMC; absID++){ // Nov 30, 2006 by PAI; was from 1 to nEMC
402       Float_t energy = 0 ;
403       // amplitude set to zero, noise will be added later
404       new((*digits)[absID]) AliEMCALDigit( -1, -1, absID, 0., TimeOfNoise(),kFALSE); // absID-1->absID
405       //look if we have to add signal?
406       
407       digit = dynamic_cast<AliEMCALDigit *>(digits->At(absID)); // absID-1->absID
408       
409       if (digit) {
410         
411         if(absID==nextSig){
412           //Add SDigits from all inputs    
413           //      ticks->Clear() ;
414           //Int_t contrib = 0 ;
415           
416           //Follow PHOS and comment out this timing model til a better one
417           //can be developed - JLK 28-Apr-2008
418           
419           //Float_t a = digit->GetAmplitude() ;
420           //Float_t b = TMath::Abs( a /fTimeSignalLength) ;
421           //Mark the beginning of the signal
422           //new((*ticks)[contrib++]) AliEMCALTick(digit->GetTime(),0, b);  
423           //Mark the end of the signal     
424           //new((*ticks)[contrib++]) AliEMCALTick(digit->GetTime()+fTimeSignalLength, -a, -b);
425           
426           // Calculate time as time of the largest digit
427           Float_t time = digit->GetTime() ;
428           Float_t aTime= digit->GetAmplitude() ;
429           
430           // loop over input
431           for(i = 0; i< fInput ; i++){  //loop over (possible) merge sources
432             if(dynamic_cast<TClonesArray *>(sdigArray->At(i))) {
433               Int_t sDigitEntries = dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast();
434               
435               if(sDigitEntries > index[i] )
436                 curSDigit = dynamic_cast<AliEMCALDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;        
437               else
438                 curSDigit = 0 ;
439               
440               //May be several digits will contribute from the same input
441               while(curSDigit && (curSDigit->GetId() == absID)){           
442                 //Shift primary to separate primaries belonging different inputs
443                 Int_t primaryoffset ;
444                 if(fManager)
445                   primaryoffset = fManager->GetMask(i) ; 
446                 else
447                   primaryoffset = i ;
448                 curSDigit->ShiftPrimary(primaryoffset) ;
449                 
450                 //Remove old timing model - JLK 28-April-2008
451                 //a = curSDigit->GetAmplitude() ;
452                 //b = a /fTimeSignalLength ;
453                 //new((*ticks)[contrib++]) AliEMCALTick(curSDigit->GetTime(),0, b);  
454                 //new((*ticks)[contrib++]) AliEMCALTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b); 
455                 if(curSDigit->GetAmplitude()>aTime) {
456                   aTime = curSDigit->GetAmplitude();
457                   time  = curSDigit->GetTime();
458                 }
459                 
460                 *digit = *digit + *curSDigit ;  //adds amplitudes of each digit
461                 
462                 index[i]++ ;
463                 if( sDigitEntries > index[i] )
464                   curSDigit = dynamic_cast<AliEMCALDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
465                 else
466                   curSDigit = 0 ;
467               }// while
468             }// source exists
469           }// loop over merging sources
470           
471           //Here we convert the summed amplitude to an energy in GeV
472           energy = sDigitizer->Calibrate(digit->GetAmplitude()) ; // GeV
473           // add fluctuations for photo-electron creation
474           energy *= static_cast<Float_t>(gRandom->Poisson(fMeanPhotonElectron)) / static_cast<Float_t>(fMeanPhotonElectron) ;
475           
476           //calculate and set time
477           //New timing model needed - JLK 28-April-2008
478           //Float_t time = FrontEdgeTime(ticks) ;
479           digit->SetTime(time) ;
480           
481           //Find next signal module
482           nextSig = nEMC + 1 ;
483           for(i = 0 ; i < fInput ; i++){
484             sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
485             if(sdigits){
486               Int_t curNext = nextSig ;
487               if(sdigits->GetEntriesFast() > index[i] && dynamic_cast<AliEMCALDigit *>(sdigits->At(index[i]))){
488                 curNext = dynamic_cast<AliEMCALDigit *>(sdigits->At(index[i]))->GetId() ;         
489               }
490               if(curNext < nextSig) nextSig = curNext ;
491             }// sdigits exist
492           } // input loop       
493           
494         }//absID==nextSig
495         
496         // add the noise now
497         energy += TMath::Abs(gRandom->Gaus(0., fPinNoise)) ;
498         // JLK 26-June-2008
499         //Now digitize the energy according to the sDigitizer method,
500         //which merely converts the energy to an integer.  Later we will
501         //check that the stored value matches our allowed dynamic ranges
502         digit->SetAmplitude(sDigitizer->Digitize(energy)) ;  
503         AliDebug(10,Form(" absID %5i energy %f nextSig %5i\n",
504                          absID, energy, nextSig));
505       }// Digit pointer exists
506       else{
507         Info("Digitizer","Digit pointer is null");
508       }
509     } // for(absID = 0; absID < nEMC; absID++)
510     
511     //ticks->Delete() ;
512     //delete ticks ;
513     
514     //JLK is it better to call Clear() here?
515     delete sdigArray ; //We should not delete its contents
516     
517     //remove digits below thresholds
518     // until 10-02-2010 remove digits with energy smaller than fDigitThreshold 3*fPinNoise
519     // now, remove digits with Digitized ADC smaller than fDigitThreshold = 3
520     Float_t energy=0;
521     for(i = 0 ; i < nEMC ; i++){
522       digit = dynamic_cast<AliEMCALDigit*>( digits->At(i) ) ;
523       if(digit){
524         //First get the energy in GeV units.
525         energy = sDigitizer->Calibrate(digit->GetAmplitude()) ;
526         //Then digitize using the calibration constants of the ocdb
527         Float_t ampADC = DigitizeEnergy(energy, digit->GetId())  ;        
528         //if(ampADC>2)printf("Digit energy %f, amp %d, amp cal %d, threshold %d\n",energy,digit->GetAmplitude(),ampADC,fDigitThreshold);
529         if(ampADC < fDigitThreshold)
530           digits->RemoveAt(i) ;
531         else 
532           digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
533       }// digit exists
534     } // digit loop
535     
536     digits->Compress() ;  
537     
538     Int_t ndigits = digits->GetEntriesFast() ; 
539     
540     //JLK 26-June-2008
541     //After we have done the summing and digitizing to create the
542     //digits, now we want to calibrate the resulting amplitude to match
543     //the dynamic range of our real data.  
544     for (i = 0 ; i < ndigits ; i++) { 
545       digit = dynamic_cast<AliEMCALDigit *>( digits->At(i) ) ; 
546       if(digit){
547         digit->SetIndexInList(i) ; 
548         energy = sDigitizer->Calibrate(digit->GetAmplitude()) ;
549         digit->SetAmplitude(DigitizeEnergy(energy, digit->GetId()) ) ;
550         //Add delay to time
551         digit->SetTime(digit->GetTime()+fTimeDelay) ;
552         // printf("digit amplitude set at end: i %d, amp %f\n",i,digit->GetAmplitude());
553       }// digit exists
554     }//Digit loop
555     
556   }//SDigitizer not null
557 }
558
559 // //_____________________________________________________________________
560 Float_t AliEMCALDigitizer::DigitizeEnergy(Float_t energy, Int_t AbsId)
561
562   // JLK 26-June-2008
563   // Returns digitized value of the energy in a cell absId
564   // using the calibration constants stored in the OCDB
565   // or default values if no CalibData object is found.
566   // This effectively converts everything to match the dynamic range
567   // of the real data we will collect
568   //
569   // Load Geometry
570   const AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
571
572   Float_t channel = -999; 
573   
574   if (geom==0){
575     AliFatal("Did not get geometry from EMCALLoader");
576   }
577   else {
578     Int_t iSupMod = -1;
579     Int_t nModule  = -1;
580     Int_t nIphi   = -1;
581     Int_t nIeta   = -1;
582     Int_t iphi    = -1;
583     Int_t ieta    = -1;
584     Bool_t bCell = geom->GetCellIndex(AbsId, iSupMod, nModule, nIphi, nIeta) ;
585     if(!bCell)
586       Error("DigitizeEnergy","Wrong cell id number : AbsId %i ", AbsId) ;
587     geom->GetCellPhiEtaIndexInSModule(iSupMod,nModule,nIphi, nIeta,iphi,ieta);
588   
589     if(fCalibData) {
590       fADCpedestalEC = fCalibData->GetADCpedestal(iSupMod,ieta,iphi);
591       fADCchannelEC = fCalibData->GetADCchannel(iSupMod,ieta,iphi);
592     }
593   
594     //channel = static_cast<Int_t> (TMath::Floor( (energy + fADCpedestalEC)/fADCchannelEC ))  ;
595     channel = (energy + fADCpedestalEC)/fADCchannelEC   ;
596
597     if(channel > fNADCEC ) 
598       channel =  fNADCEC ; 
599   }// geometry exists
600   
601   return channel ;
602   
603 }
604
605 //____________________________________________________________________________
606 void AliEMCALDigitizer::Exec(Option_t *option) 
607
608   // Steering method to process digitization for events
609   // in the range from fFirstEvent to fLastEvent.
610   // This range is optionally set by SetEventRange().
611   // if fLastEvent=-1, then process events until the end.
612   // by default fLastEvent = fFirstEvent (process only one event)
613
614   if (!fInit) { // to prevent overwrite existing file
615     Error( "Exec", "Give a version name different from %s", fEventFolderName.Data() ) ;
616     return ;
617   } 
618
619   if (strstr(option,"print")) {
620
621     Print();
622     return ; 
623   }
624   
625   if(strstr(option,"tim"))
626     gBenchmark->Start("EMCALDigitizer");
627
628   AliRunLoader *rl = AliRunLoader::Instance();
629   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
630   Int_t nEvents = 0;
631   if(!emcalLoader){
632     AliFatal("Did not get the  Loader");
633   }
634   else{
635     // Post Digitizer to the white board
636     emcalLoader->PostDigitizer(this) ;
637     
638     if (fLastEvent == -1)  {
639       fLastEvent = rl->GetNumberOfEvents() - 1 ;
640     }
641     else if (fManager) 
642       fLastEvent = fFirstEvent ; // what is this ??
643     
644     nEvents = fLastEvent - fFirstEvent + 1;
645     Int_t ievent  = -1;
646
647     TClonesArray* digitsTMP = new TClonesArray("AliEMCALDigit",    32*96);
648     TClonesArray* digitsTRG = new TClonesArray("AliEMCALRawDigit", 32*96);
649
650     rl->LoadSDigits("EMCAL");
651     for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
652       
653       rl->GetEvent(ievent);
654       
655       Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
656       
657       WriteDigits() ;
658       
659       //Trigger Digits
660       //-------------------------------------
661      
662
663       Digits2FastOR(digitsTMP, digitsTRG);  
664       
665       WriteDigits(digitsTRG);
666       
667       (emcalLoader->TreeD())->Fill();
668       
669       emcalLoader->WriteDigits(   "OVERWRITE");
670       emcalLoader->WriteDigitizer("OVERWRITE");
671       
672       Unload();
673       
674       digitsTRG  ->Delete();
675       digitsTMP  ->Delete();
676       
677       //-------------------------------------
678       
679       if(strstr(option,"deb"))
680         PrintDigits(option);
681       if(strstr(option,"table")) gObjectTable->Print();
682       
683       //increment the total number of Digits per run 
684       fDigitsInRun += emcalLoader->Digits()->GetEntriesFast() ;  
685     }//loop
686     
687     emcalLoader->CleanDigitizer() ;
688     
689   }//loader exists
690   
691   if(strstr(option,"tim")){
692     gBenchmark->Stop("EMCALDigitizer");
693     AliInfo(Form("Exec: took %f seconds for Digitizing %f seconds per event", 
694          gBenchmark->GetCpuTime("EMCALDigitizer"), gBenchmark->GetCpuTime("EMCALDigitizer")/nEvents )) ;
695   } 
696 }
697
698 //____________________________________________________________________________ 
699 void AliEMCALDigitizer::Digits2FastOR(TClonesArray* digitsTMP, TClonesArray* digitsTRG)
700 {
701   // FEE digits afterburner to produce TRG digits 
702   // we are only interested in the FEE digit deposited energy
703   // to be converted later into a voltage value
704   
705   // push the FEE digit to its associated FastOR (numbered from 0:95)
706   // TRU is in charge of summing module digits
707   
708   AliRunLoader *runLoader = AliRunLoader::Instance();
709   
710   AliRun* run = runLoader->GetAliRun();
711   
712   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(runLoader->GetDetectorLoader("EMCAL"));
713   if(!emcalLoader){
714     AliFatal("Did not get the  Loader");
715   }
716   else {
717     AliEMCALGeometry* geom = dynamic_cast<AliEMCAL*>(run->GetDetector("EMCAL"))->GetGeometry();
718     
719     // build FOR from simulated digits
720     // and xfer to the corresponding TRU input (mapping)
721     
722     TClonesArray* digits = emcalLoader->Digits();
723     
724     TIter NextDigit(digits);
725     while (AliEMCALDigit* digit = (AliEMCALDigit*)NextDigit())
726     {
727       Int_t id = digit->GetId();
728       
729       Int_t trgid;
730       if (geom && geom->GetFastORIndexFromCellIndex(id , trgid)) 
731       {
732         AliDebug(1,Form("trigger digit id: %d from cell id: %d\n",trgid,id));
733         
734         AliEMCALDigit* d = static_cast<AliEMCALDigit*>(digitsTMP->At(trgid));
735         
736         if (!d)
737               {
738           new((*digitsTMP)[trgid]) AliEMCALDigit(*digit);
739           d = (AliEMCALDigit*)digitsTMP->At(trgid);
740           d->SetId(trgid);
741               } 
742         else
743               {
744           *d = *d + *digit;
745               }
746       }
747     }
748     
749     if (AliDebugLevel()) printf("Number of TRG digits: %d\n",digitsTMP->GetEntriesFast());
750     
751     Int_t    nSamples = 32;
752     Int_t *timeSamples = new Int_t[nSamples];
753     
754     NextDigit = TIter(digitsTMP);
755     while (AliEMCALDigit* digit = (AliEMCALDigit*)NextDigit())
756     {
757       if (digit)
758       {
759         Int_t     id = digit->GetId();
760         Float_t time = digit->GetTime();
761         
762         Double_t depositedEnergy = 0.;
763         for (Int_t j = 1; j <= digit->GetNprimary(); j++) depositedEnergy += digit->GetDEPrimary(j);
764         
765         if (AliDebugLevel()) printf("Deposited Energy: %f\n", depositedEnergy);
766         
767         // FIXME: Check digit time!
768         if (depositedEnergy)
769               {
770           DigitalFastOR(time, depositedEnergy, timeSamples, nSamples);
771           
772           for (Int_t j=0;j<nSamples;j++) 
773           {
774             timeSamples[j] = ((j << 12) & 0xFF000) | (timeSamples[j] & 0xFFF);
775           }
776           
777           new((*digitsTRG)[digitsTRG->GetEntriesFast()]) AliEMCALRawDigit(id, timeSamples, nSamples);
778               }
779       }
780     }
781     
782     delete [] timeSamples;
783     
784   }// loader exists
785   
786 }
787
788 //____________________________________________________________________________ 
789 void AliEMCALDigitizer::DigitalFastOR( Double_t time, Double_t dE, Int_t timeSamples[], Int_t nSamples )
790 {
791         // parameters:  
792         // id: 0..95
793         const Int_t    reso = 11;      // 11-bit resolution ADC
794         const Double_t vFSR = 1;       // Full scale input voltage range
795         const Double_t Ne   = 125;     // signal of the APD per MeV of energy deposit in a tower: 125 photo-e-/MeV @ M=30
796         const Double_t vA   = .136e-6; // CSP output range: 0.136uV/e-
797         const Double_t rise = 40e-9;   // rise time (10-90%) of the FastOR signal before shaping
798         
799         const Double_t kTimeBinWidth = 25E-9; // sampling frequency (40MHz)
800         
801         Double_t vV = 1000. * dE * Ne * vA; // GeV 2 MeV
802         
803         TF1 signalF("signal", AnalogFastORFunction, 0, nSamples * kTimeBinWidth, 3);
804         signalF.SetParameter( 0,   vV ); 
805         signalF.SetParameter( 1, time ); // FIXME: when does the signal arrive? Might account for cable lengths
806         signalF.SetParameter( 2, rise );
807         
808         for (Int_t iTime=0; iTime<nSamples; iTime++) 
809         {
810                 // FIXME: add noise (probably not simply Gaussian) according to DA measurements
811                 // probably plan an access to OCDB
812                 
813                 timeSamples[iTime] = int((TMath::Power(2, reso) / vFSR) * signalF.Eval(iTime * kTimeBinWidth) + 0.5);
814         }
815 }
816
817
818 //____________________________________________________________________________ 
819 //Float_t AliEMCALDigitizer::FrontEdgeTime(TClonesArray * ticks) 
820 //{ 
821 //  //  Returns the shortest time among all time ticks
822 //
823 //  ticks->Sort() ; //Sort in accordance with times of ticks
824 //  TIter it(ticks) ;
825 //  AliEMCALTick * ctick = (AliEMCALTick *) it.Next() ;
826 //  Float_t time = ctick->CrossingTime(fTimeThreshold) ;    
827 //  
828 //  AliEMCALTick * t ;  
829 //  while((t=(AliEMCALTick*) it.Next())){
830 //    if(t->GetTime() < time)  //This tick starts before crossing
831 //      *ctick+=*t ;
832 //    else
833 //      return time ;
834 //    
835 //    time = ctick->CrossingTime(fTimeThreshold) ;    
836 //  }
837 //  return time ;
838 //}
839 //
840
841 //____________________________________________________________________________ 
842 Bool_t AliEMCALDigitizer::Init()
843 {
844   // Makes all memory allocations
845   fInit = kTRUE ; 
846   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
847
848   if ( emcalLoader == 0 ) {
849     Fatal("Init", "Could not obtain the AliEMCALLoader");  
850     return kFALSE;
851   } 
852
853   fFirstEvent = 0 ; 
854   fLastEvent = fFirstEvent ; 
855   
856   if(fManager)
857     fInput = fManager->GetNinputs() ; 
858   else 
859     fInput           = 1 ;
860
861   fInputFileNames  = new TString[fInput] ; 
862   fEventNames      = new TString[fInput] ; 
863   fInputFileNames[0] = GetTitle() ; 
864   fEventNames[0]     = fEventFolderName.Data() ; 
865   Int_t index ; 
866   for (index = 1 ; index < fInput ; index++) {
867     fInputFileNames[index] = dynamic_cast<AliStream*>(fManager->GetInputStream(index))->GetFileName(0); 
868     TString tempo = fManager->GetInputFolderName(index) ;
869     fEventNames[index] = tempo.Remove(tempo.Length()-1) ; // strip of the stream number added bt fManager 
870   }
871   
872   //to prevent cleaning of this object while GetEvent is called
873   emcalLoader->GetDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
874
875   //Calibration instance
876   fCalibData = emcalLoader->CalibData();
877   return fInit ;    
878 }
879
880 //____________________________________________________________________________ 
881 void AliEMCALDigitizer::InitParameters()
882
883   // Parameter initialization for digitizer
884   
885   // Get the parameters from the OCDB via the loader
886   AliRunLoader *rl = AliRunLoader::Instance();
887   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
888   AliEMCALSimParam * simParam = 0x0;
889   if(emcalLoader) simParam = emcalLoader->SimulationParameters();
890         
891   if(!simParam){
892           simParam = AliEMCALSimParam::GetInstance();
893           AliWarning("Simulation Parameters not available in OCDB?");
894   }
895         
896   fMeanPhotonElectron = simParam->GetMeanPhotonElectron();//4400;  // electrons per GeV 
897   fPinNoise           = simParam->GetPinNoise();//0.012; // pin noise in GeV from analysis test beam data 
898   if (fPinNoise < 0.0001 ) 
899     Warning("InitParameters", "No noise added\n") ; 
900   fDigitThreshold     = simParam->GetDigitThreshold(); //fPinNoise * 3; // 3 * sigma
901   fTimeResolution     = simParam->GetTimeResolution(); //0.6e-9 ; // 600 pc
902   fTimeDelay          = simParam->GetTimeDelay(); //600e-9 ; // 600 ns
903
904   // These defaults are normally not used. 
905   // Values are read from calibration database instead
906   fADCchannelEC       = 0.0153; // Update 24 Apr 2007: 250./16/1024 - width of one ADC channel in GeV
907   fADCpedestalEC      = 0.0 ;  // GeV
908
909   fNADCEC          = simParam->GetNADCEC();//(Int_t) TMath::Power(2,16) ;  // number of channels in Tower ADC - 65536
910
911   AliDebug(2,Form("Mean Photon Electron %d, Noise %f, Digit Threshold %d,Time Resolution %g,NADCEC %d",
912                 fMeanPhotonElectron,fPinNoise,fDigitThreshold,fTimeResolution,fNADCEC));
913
914   // Not used anymore, remove?
915   // fTimeSignalLength   = 1.0e-9 ;
916   // fTimeThreshold      = 0.001*10000000 ; // Means 1 MeV in terms of SDigits amplitude ??
917
918 }
919
920 //__________________________________________________________________
921 void AliEMCALDigitizer::MixWith(TString alirunFileName, TString eventFolderName)
922 {
923   // Allows to produce digits by superimposing background and signal event.
924   // It is assumed, that headers file with SIGNAL events is opened in 
925   // the constructor. 
926   // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed 
927   // Thus we avoid writing (changing) huge and expensive 
928   // backgound files: all output will be writen into SIGNAL, i.e. 
929   // opened in constructor file. 
930   //
931   // One can open as many files to mix with as one needs.
932   // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
933   // can be mixed.
934
935   if( strcmp(GetName(), "") == 0 )
936     Init() ;
937   
938   if(fManager){
939     Error("MixWith", "Cannot use this method under AliRunDigitizer") ;
940     return ;
941   } 
942   // looking for file which contains AliRun
943   if (gSystem->AccessPathName(alirunFileName)) {// file does not exist
944     Error("MixWith", "File %s does not exist!", alirunFileName.Data()) ;
945     return ; 
946   }
947   // looking for the file which contains SDigits
948   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
949   
950   if(!emcalLoader){
951     AliFatal("Did not get the  Loader");
952   }
953   else{
954     TString fileName( emcalLoader->GetSDigitsFileName() ) ; 
955     if ( eventFolderName != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name 
956       fileName = fileName.ReplaceAll(".root", "") + "_" + eventFolderName + ".root" ;
957     if ( (gSystem->AccessPathName(fileName)) ) { 
958       Error("MixWith", "The file %s does not exist!", fileName.Data()) ;
959       return ;
960     }
961     // need to increase the arrays
962     // MvL: This code only works when fInput == 1, I think.
963     TString tempo = fInputFileNames[fInput-1] ; 
964     delete [] fInputFileNames ; 
965     fInputFileNames = new TString[fInput+1] ; 
966     fInputFileNames[fInput-1] = tempo ; 
967     
968     tempo = fEventNames[fInput-1] ; 
969     delete [] fEventNames ; 
970     fEventNames = new TString[fInput+1] ; 
971     fEventNames[fInput-1] = tempo ; 
972     
973     fInputFileNames[fInput] = alirunFileName ; 
974     fEventNames[fInput]     = eventFolderName ;
975     fInput++ ;
976     
977   }// Loader exits
978   
979 }  
980
981 //__________________________________________________________________
982 void AliEMCALDigitizer::Print1(Option_t * option)
983 { // 19-nov-04 - just for convinience
984   Print(); 
985   PrintDigits(option);
986 }
987
988 //__________________________________________________________________
989 void AliEMCALDigitizer::Print(Option_t*)const 
990 {
991   // Print Digitizer's parameters
992   printf("Print: \n------------------- %s -------------", GetName() ) ; 
993   if( strcmp(fEventFolderName.Data(), "") != 0 ){
994     printf(" Writing Digits to branch with title  %s\n", fEventFolderName.Data()) ;
995     
996     Int_t nStreams ; 
997     if (fManager) 
998       nStreams =  GetNInputStreams() ;
999     else 
1000       nStreams = fInput ; 
1001   
1002     AliRunLoader *rl=0;
1003     
1004     Int_t index = 0 ;  
1005     for (index = 0 ; index < nStreams ; index++) {  
1006       TString tempo(fEventNames[index]) ; 
1007       tempo += index ;
1008       if ((rl = AliRunLoader::GetRunLoader(tempo)) == 0)
1009         rl = AliRunLoader::Open(fInputFileNames[index], tempo) ; 
1010       AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
1011       if(emcalLoader){
1012         TString fileName( emcalLoader->GetSDigitsFileName() ) ; 
1013         if ( fEventNames[index] != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name 
1014           fileName = fileName.ReplaceAll(".root", "") + "_" + fEventNames[index]  + ".root" ;
1015         printf ("Adding SDigits from %s %s\n", fInputFileNames[index].Data(), fileName.Data()) ; 
1016       }//loader
1017     }
1018
1019     AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1020
1021     if(emcalLoader) printf("\nWriting digits to %s", emcalLoader->GetDigitsFileName().Data()) ;
1022     else printf("\nNULL LOADER");
1023
1024     printf("\nWith following parameters:\n") ;
1025     
1026     printf("    Electronics noise in EMC (fPinNoise) = %f\n", fPinNoise) ;
1027     printf("    Threshold  in Tower  (fDigitThreshold) = %d\n", fDigitThreshold)  ;
1028     printf("---------------------------------------------------\n")  ;
1029   }
1030   else
1031     printf("Print: AliEMCALDigitizer not initialized") ; 
1032 }
1033
1034 //__________________________________________________________________
1035 void AliEMCALDigitizer::PrintDigits(Option_t * option)
1036 {
1037   //utility method for printing digit information
1038
1039   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1040   if(emcalLoader){
1041     TClonesArray * digits  = emcalLoader->Digits() ;
1042     TClonesArray * sdigits = emcalLoader->SDigits() ;
1043     
1044     printf("\n #Digits: %d : sdigits %d ", digits->GetEntriesFast(), sdigits->GetEntriesFast()) ; 
1045     printf("\n event %d", emcalLoader->GetRunLoader()->GetEventNumber());
1046     
1047     if(strstr(option,"all")){  
1048       //loop over digits
1049       AliEMCALDigit * digit;
1050       printf("\nEMC digits (with primaries):\n")  ;
1051       printf("\n   Id  Amplitude    Time          Index Nprim: Primaries list \n") ;    
1052       Int_t index ;
1053       for (index = 0 ; index < digits->GetEntries() ; index++) {
1054         digit = dynamic_cast<AliEMCALDigit *>(digits->At(index)) ;
1055         if(digit){
1056           printf("\n%6d  %8f    %6.5e %4d      %2d : ",
1057                  digit->GetId(), digit->GetAmplitude(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;  
1058           Int_t iprimary;
1059           for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
1060             printf("%d ",digit->GetPrimary(iprimary+1) ) ; 
1061           }
1062         }// digit exists
1063       }// loop
1064     }
1065     printf("\n");
1066   }// loader exists
1067   else printf("NULL LOADER, cannot print\n");
1068 }
1069
1070 //__________________________________________________________________
1071 Float_t AliEMCALDigitizer::TimeOfNoise(void)
1072 {  
1073   // Calculates the time signal generated by noise
1074   //PH  Info("TimeOfNoise", "Change me") ; 
1075   return gRandom->Rndm() * 1.28E-5;
1076 }
1077
1078 //__________________________________________________________________
1079 void AliEMCALDigitizer::Unload() 
1080 {  
1081   // Unloads the SDigits and Digits
1082   AliRunLoader *rl=0;
1083     
1084   Int_t i ; 
1085   for(i = 1 ; i < fInput ; i++){
1086     TString tempo(fEventNames[i]) ; 
1087     tempo += i;
1088     if ((rl = AliRunLoader::GetRunLoader(tempo))) 
1089       rl->GetDetectorLoader("EMCAL")->UnloadSDigits() ; 
1090   }
1091   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1092   if(emcalLoader)emcalLoader->UnloadDigits() ; 
1093   else AliFatal("Did not get the loader");
1094 }
1095
1096 //_________________________________________________________________________________________
1097 void AliEMCALDigitizer::WriteDigits()
1098 {
1099   
1100   // Makes TreeD in the output file. 
1101   // Check if branch already exists: 
1102   //   if yes, exit without writing: ROOT TTree does not support overwriting/updating of 
1103   //      already existing branches. 
1104   //   else creates branch with Digits, named "EMCAL", title "...",
1105   //      and branch "AliEMCALDigitizer", with the same title to keep all the parameters
1106   //      and names of files, from which digits are made.
1107   
1108   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1109   
1110   if(emcalLoader){
1111     const TClonesArray * digits = emcalLoader->Digits() ; 
1112     TTree * treeD = emcalLoader->TreeD(); 
1113     if ( !treeD ) {
1114       emcalLoader->MakeDigitsContainer();
1115       treeD = emcalLoader->TreeD(); 
1116     }
1117     
1118     // -- create Digits branch
1119     Int_t bufferSize = 32000 ;    
1120     TBranch * digitsBranch = 0;
1121     if ((digitsBranch = treeD->GetBranch("EMCAL"))) {
1122       digitsBranch->SetAddress(&digits);
1123       AliWarning("Digits Branch already exists. Not all digits will be visible");
1124     }
1125     else
1126       treeD->Branch("EMCAL","TClonesArray",&digits,bufferSize);
1127     //digitsBranch->SetTitle(fEventFolderName);
1128     
1129     //  treeD->Fill() ;
1130     /*  
1131      emcalLoader->WriteDigits("OVERWRITE");
1132      emcalLoader->WriteDigitizer("OVERWRITE");
1133      
1134      Unload() ; 
1135      */
1136     
1137   }// loader exists
1138   else AliFatal("Loader not available");
1139 }
1140
1141 //__________________________________________________________________
1142 void AliEMCALDigitizer::WriteDigits(TClonesArray* digits, const char* branchName)
1143 {
1144         //
1145         AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::Instance()->GetDetectorLoader("EMCAL"));
1146         if(emcalLoader){
1147     
1148     TTree* treeD = emcalLoader->TreeD(); 
1149     if (!treeD) 
1150     {
1151       emcalLoader->MakeDigitsContainer();
1152       treeD = emcalLoader->TreeD(); 
1153     }
1154     
1155     // -- create Digits branch
1156     Int_t bufferSize = 32000;
1157     
1158     if (TBranch* triggerBranch = treeD->GetBranch(branchName)) 
1159     {
1160       triggerBranch->SetAddress(&digits);
1161     }
1162     else
1163     {
1164       treeD->Branch(branchName,"TClonesArray",&digits,bufferSize);
1165     }
1166     
1167     //  treeD->Fill();
1168   }// loader exists
1169   else AliFatal("Loader not available");
1170 }
1171
1172 //__________________________________________________________________
1173 void AliEMCALDigitizer::Browse(TBrowser* b)
1174 {
1175   TTask::Browse(b);
1176 }