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