]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALDigitizer.cxx
Restoring EMCAL digitization and reconstruction
[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 
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 "TList.h"
67 #include "TH1.h"
68 #include "TBrowser.h"
69 #include "TObjectTable.h"
70
71 // --- AliRoot header files ---
72 #include "AliRunDigitizer.h"
73 #include "AliEMCALDigit.h"
74 #include "AliEMCAL.h"
75 #include "AliEMCALGetter.h"
76 #include "AliEMCALDigitizer.h"
77 #include "AliEMCALSDigitizer.h"
78 #include "AliEMCALGeometry.h"
79 #include "AliEMCALTick.h"
80 #include "AliEMCALJetMicroDst.h"
81
82 ClassImp(AliEMCALDigitizer)
83
84
85 //____________________________________________________________________________ 
86   AliEMCALDigitizer::AliEMCALDigitizer():AliDigitizer("",""),
87                                        fInput(0),
88                                        fInputFileNames(0x0),
89                                        fEventNames(0x0)
90 {
91   // ctor
92   InitParameters() ; 
93   fDefaultInit = kTRUE ; 
94   fManager = 0 ;                     // We work in the standalong mode
95   fEventFolderName = "" ; 
96 }
97
98 //____________________________________________________________________________ 
99 AliEMCALDigitizer::AliEMCALDigitizer(TString alirunFileName, TString eventFolderName):
100   AliDigitizer("EMCAL"+AliConfig::Instance()->GetDigitizerTaskName(), alirunFileName),
101   fInputFileNames(0), fEventNames(0), fEventFolderName(eventFolderName)
102 {
103   // ctor
104
105   InitParameters() ; 
106   Init() ;
107   fDefaultInit = kFALSE ; 
108   fManager = 0 ;                     // We work in the standalong mode
109 }
110
111 //____________________________________________________________________________ 
112 AliEMCALDigitizer::AliEMCALDigitizer(const AliEMCALDigitizer & d) : AliDigitizer(d)
113 {
114   // copyy ctor 
115
116   SetName(d.GetName()) ; 
117   SetTitle(d.GetTitle()) ; 
118   fDigitThreshold = d.fDigitThreshold ; 
119   fMeanPhotonElectron = d.fMeanPhotonElectron ; 
120   fPedestal           = d.fPedestal ; 
121   fSlope              = d.fSlope ; 
122   fPinNoise           = d.fPinNoise ; 
123   fTimeResolution     = d.fTimeResolution ; 
124   fTimeThreshold      = d.fTimeThreshold ; 
125   fTimeSignalLength   = d.fTimeSignalLength ; 
126   fADCchannelEC       = d.fADCchannelEC ; 
127   fADCpedestalEC      = d.fADCpedestalEC ; 
128   fNADCEC             = d.fNADCEC ;
129   fEventFolderName    = d.fEventFolderName;
130  }
131
132 //____________________________________________________________________________ 
133 AliEMCALDigitizer::AliEMCALDigitizer(AliRunDigitizer * rd):
134  AliDigitizer(rd,"EMCAL"+AliConfig::Instance()->GetDigitizerTaskName()),
135  fEventFolderName(0)
136 {
137   // ctor Init() is called by RunDigitizer
138   fManager = rd ; 
139   fEventFolderName = fManager->GetInputFolderName(0) ;
140   SetTitle(dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetFileName(0));
141   InitParameters() ; 
142   fDefaultInit = kFALSE ; 
143 }
144
145 //____________________________________________________________________________ 
146   AliEMCALDigitizer::~AliEMCALDigitizer()
147 {
148   AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle()) ;
149    
150   // Clean Digitizer from the white board
151   gime->EmcalLoader()->CleanDigitizer() ;
152
153   delete [] fInputFileNames ; 
154   delete [] fEventNames ; 
155
156 }
157
158 //____________________________________________________________________________
159 void AliEMCALDigitizer::Digitize(Int_t event) 
160
161
162   // Makes the digitization of the collected summable digits
163   // for this it first creates the array of all EMCAL modules
164   // filled with noise (different for EMC, CPV and PPSD) and
165   // after that adds contributions from SDigits. This design 
166   // helps to avoid scanning over the list of digits to add 
167   // contribution of any new SDigit.
168   static int isTrd1Geom = -1; // -1 - mean undefined 
169   static int nEMC=0; //max number of digits possible
170
171   AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle()) ; 
172   Int_t ReadEvent = event ; 
173   // fManager is data member from AliDigitizer
174   if (fManager) 
175     ReadEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(0))->GetCurrentEventNumber() ; 
176   AliDebug(1,Form("Adding event %d from input stream 0 %s %s", 
177                   ReadEvent, GetTitle(), fEventFolderName.Data())) ; 
178   gime->Event(ReadEvent, "S") ;
179   TClonesArray * digits = gime->Digits() ; 
180   digits->Clear() ;
181
182   const AliEMCALGeometry *geom = gime->EMCALGeometry() ;
183   if(isTrd1Geom < 0) { 
184     TString ng(geom->GetName());
185     isTrd1Geom = 0;
186     if(ng.Contains("SHISH") &&  ng.Contains("TRD1")) isTrd1Geom = 1;
187
188     if(isTrd1Geom == 0) nEMC = geom->GetNPhi()*geom->GetNZ();
189     else                nEMC = geom->GetNCells();
190     printf(" nEMC %i (number cells in EMCAL) | %s | isTrd1Geom %i\n", nEMC, geom->GetName(), isTrd1Geom);
191   }
192   Int_t absID ;
193
194   digits->Expand(nEMC) ;
195
196   // get first the sdigitizer from the tasks list (must have same name as the digitizer)
197   if ( !gime->SDigitizer() ) 
198     gime->LoadSDigitizer();
199   AliEMCALSDigitizer * sDigitizer = gime->SDigitizer(); 
200   
201   if ( !sDigitizer )
202     Fatal("Digitize", "SDigitizer with name %s %s not found", fEventFolderName.Data(), GetTitle() ) ; 
203
204   //take all the inputs to add together and load the SDigits
205   TObjArray * sdigArray = new TObjArray(fInput) ;
206   sdigArray->AddAt(gime->SDigits(), 0) ;
207   Int_t i ;
208   for(i = 1 ; i < fInput ; i++){
209     TString tempo(fEventNames[i]) ; 
210     tempo += i ;
211     AliEMCALGetter * gime = AliEMCALGetter::Instance(fInputFileNames[i], tempo) ; 
212     if (fManager) 
213       ReadEvent = dynamic_cast<AliStream*>(fManager->GetInputStream(i))->GetCurrentEventNumber() ; 
214     Info("Digitize", "Adding event %d from input stream %d %s %s", ReadEvent, i, fInputFileNames[i].Data(), tempo.Data()) ; 
215     gime->Event(ReadEvent,"S");
216     sdigArray->AddAt(gime->SDigits(), i) ;
217   }
218  
219   //Find the first tower with signal
220   Int_t nextSig = nEMC + 1 ; 
221   TClonesArray * sdigits ;  
222   for(i = 0 ; i < fInput ; i++){
223     sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
224     if ( !sdigits->GetEntriesFast() )
225       continue ; 
226     Int_t curNext = dynamic_cast<AliEMCALDigit *>(sdigits->At(0))->GetId() ;
227      if(curNext < nextSig) 
228        nextSig = curNext ;
229      AliDebug(1,Form("input %i : #sdigits %i \n",
230                      i, sdigits->GetEntriesFast()));
231   }
232   AliDebug(1,Form("FIRST tower with signal %i \n", nextSig));
233
234   TArrayI index(fInput) ;
235   index.Reset() ;  //Set all indexes to zero
236
237   AliEMCALDigit * digit ;
238   AliEMCALDigit * curSDigit ;
239
240   TClonesArray * ticks = new TClonesArray("AliEMCALTick",1000) ;
241
242   //Put Noise contributiony
243   for(absID = 1; absID <= nEMC; absID++){
244     Float_t amp = 0 ;
245     // amplitude set to zero, noise will be added later
246     new((*digits)[absID-1]) AliEMCALDigit( -1, -1, absID, 0, TimeOfNoise() ) ;
247     //look if we have to add signal?
248     digit = dynamic_cast<AliEMCALDigit *>(digits->At(absID-1)) ;
249     
250     if(absID==nextSig){
251       //Add SDigits from all inputs    
252       ticks->Clear() ;
253       Int_t contrib = 0 ;
254       Float_t a = digit->GetAmp() ;
255       Float_t b = TMath::Abs( a /fTimeSignalLength) ;
256       //Mark the beginning of the signal
257       new((*ticks)[contrib++]) AliEMCALTick(digit->GetTime(),0, b);  
258       //Mark the end of the signal     
259       new((*ticks)[contrib++]) AliEMCALTick(digit->GetTime()+fTimeSignalLength, -a, -b);
260       
261       // loop over input
262       for(i = 0; i< fInput ; i++){  //loop over (possible) merge sources
263         if(dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
264           curSDigit = dynamic_cast<AliEMCALDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;      
265         else
266           curSDigit = 0 ;
267         //May be several digits will contribute from the same input
268         while(curSDigit && (curSDigit->GetId() == absID)){         
269           //Shift primary to separate primaries belonging different inputs
270           Int_t primaryoffset ;
271           if(fManager)
272             primaryoffset = fManager->GetMask(i) ; 
273           else
274             primaryoffset = i ;
275           curSDigit->ShiftPrimary(primaryoffset) ;
276           
277           a = curSDigit->GetAmp() ;
278           b = a /fTimeSignalLength ;
279           new((*ticks)[contrib++]) AliEMCALTick(curSDigit->GetTime(),0, b);  
280           new((*ticks)[contrib++]) AliEMCALTick(curSDigit->GetTime()+fTimeSignalLength, -a, -b); 
281
282           *digit = *digit + *curSDigit ;  //add energies
283
284           index[i]++ ;
285           if( dynamic_cast<TClonesArray *>(sdigArray->At(i))->GetEntriesFast() > index[i] )
286             curSDigit = dynamic_cast<AliEMCALDigit*>(dynamic_cast<TClonesArray *>(sdigArray->At(i))->At(index[i])) ;
287           else
288             curSDigit = 0 ;
289         }
290       }
291       // add fluctuations for photo-electron creation
292       amp = sDigitizer->Calibrate(digit->GetAmp()) ; // GeV
293       amp *= static_cast<Float_t>(gRandom->Poisson(fMeanPhotonElectron)) / static_cast<Float_t>(fMeanPhotonElectron) ;
294   
295       //calculate and set time
296       Float_t time = FrontEdgeTime(ticks) ;
297       digit->SetTime(time) ;
298
299       //Find next signal module
300       nextSig = nEMC + 1 ;
301       for(i = 0 ; i < fInput ; i++){
302         sdigits = dynamic_cast<TClonesArray *>(sdigArray->At(i)) ;
303         Int_t curNext = nextSig ;
304         if(sdigits->GetEntriesFast() > index[i] ){
305           curNext = dynamic_cast<AliEMCALDigit *>(sdigits->At(index[i]))->GetId() ;       
306         }
307         if(curNext < nextSig) nextSig = curNext ;
308       }
309     }
310     // add the noise now
311     amp += TMath::Abs(gRandom->Gaus(0., fPinNoise)) ;
312     digit->SetAmp(sDigitizer->Digitize(amp)) ;  
313     AliDebug(10,Form(" absID %5i amp %f nextSig %5i\n",
314                      absID, amp, nextSig));
315   } // for(absID = 1; absID <= nEMC; absID++)
316   
317   ticks->Delete() ;
318   delete ticks ;
319
320   delete sdigArray ; //We should not delete its contents
321
322   //remove digits below thresholds
323   for(i = 0 ; i < nEMC ; i++){
324     digit = dynamic_cast<AliEMCALDigit*>( digits->At(i) ) ;
325     Float_t threshold = fDigitThreshold ; 
326     if(sDigitizer->Calibrate( digit->GetAmp() ) < threshold)
327       digits->RemoveAt(i) ;
328     else 
329       digit->SetTime(gRandom->Gaus(digit->GetTime(),fTimeResolution) ) ;
330   }
331   
332   digits->Compress() ;  
333   
334   Int_t ndigits = digits->GetEntriesFast() ; 
335   digits->Expand(ndigits) ;
336   
337   //Set indexes in list of digits and fill hists.
338   sv::FillH1(fHists, 0, Double_t(ndigits));
339   Float_t energy=0., esum=0.;
340   for (i = 0 ; i < ndigits ; i++) { 
341     digit = dynamic_cast<AliEMCALDigit *>( digits->At(i) ) ; 
342     digit->SetIndexInList(i) ; 
343     energy = sDigitizer->Calibrate(digit->GetAmp()) ;
344     esum += energy;
345     digit->SetAmp(DigitizeEnergy(energy) ) ; // for what ??
346     sv::FillH1(fHists, 2, double(digit->GetAmp()));
347     sv::FillH1(fHists, 3, double(energy));
348     sv::FillH1(fHists, 4, double(digit->GetId()));
349   }
350   sv::FillH1(fHists, 1, esum);
351 }
352
353 //____________________________________________________________________________
354
355 Int_t AliEMCALDigitizer::DigitizeEnergy(Float_t energy)
356
357   Int_t channel = -999; 
358   channel = static_cast<Int_t> (TMath::Ceil( (energy + fADCpedestalEC)/fADCchannelEC ))  ;
359   if(channel > fNADCEC ) 
360     channel =  fNADCEC ; 
361   return channel ;
362 }
363
364 //____________________________________________________________________________
365 void AliEMCALDigitizer::Exec(Option_t *option) 
366
367   // Steering method to process digitization for events
368   // in the range from fFirstEvent to fLastEvent.
369   // This range is optionally set by SetEventRange().
370   // if fLastEvent=-1, then process events until the end.
371   // by default fLastEvent = fFirstEvent (process only one event)
372
373   if (!fInit) { // to prevent overwrite existing file
374     Error( "Exec", "Give a version name different from %s", fEventFolderName.Data() ) ;
375     return ;
376   } 
377
378   if (strstr(option,"print")) {
379     Print();
380     return ; 
381   }
382   
383   if(strstr(option,"tim"))
384     gBenchmark->Start("EMCALDigitizer");
385   
386   AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle()) ;
387    
388   // Post Digitizer to the white board
389   gime->PostDigitizer(this) ;
390   
391   if (fLastEvent == -1) 
392     fLastEvent = gime->MaxEvent() - 1 ;
393   else if (fManager) 
394     fLastEvent = fFirstEvent ; // what is this ??
395
396   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
397   Int_t ievent;
398
399   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
400     gime->Event(ievent,"S") ; 
401
402     Digitize(ievent) ; //Add prepared SDigits to digits and add the noise
403
404     WriteDigits() ;
405
406     if(strstr(option,"deb"))
407       PrintDigits(option);
408     if(strstr(option,"table")) gObjectTable->Print();
409
410     //increment the total number of Digits per run 
411     fDigitsInRun += gime->Digits()->GetEntriesFast() ;  
412   }
413   //  gime->WriteDigitizer("OVERWRITE");
414
415   gime->EmcalLoader()->CleanDigitizer() ;
416
417   if(strstr(option,"tim")){
418     gBenchmark->Stop("EMCALDigitizer");
419     printf("Exec: took %f seconds for Digitizing %f seconds per event", 
420          gBenchmark->GetCpuTime("EMCALDigitizer"), gBenchmark->GetCpuTime("EMCALDigitizer")/nEvents ) ;
421   } 
422 }
423
424 //____________________________________________________________________________ 
425 Float_t AliEMCALDigitizer::FrontEdgeTime(TClonesArray * ticks) 
426
427   //  Returns the shortest time among all time ticks
428
429   ticks->Sort() ; //Sort in accordance with times of ticks
430   TIter it(ticks) ;
431   AliEMCALTick * ctick = (AliEMCALTick *) it.Next() ;
432   Float_t time = ctick->CrossingTime(fTimeThreshold) ;    
433   
434   AliEMCALTick * t ;  
435   while((t=(AliEMCALTick*) it.Next())){
436     if(t->GetTime() < time)  //This tick starts before crossing
437       *ctick+=*t ;
438     else
439       return time ;
440     
441     time = ctick->CrossingTime(fTimeThreshold) ;    
442   }
443   return time ;
444 }
445
446 //____________________________________________________________________________ 
447 Bool_t AliEMCALDigitizer::Init()
448 {
449   // Makes all memory allocations
450   fInit = kTRUE ; 
451   AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName) ;  
452   if ( gime == 0 ) {
453     Error("Init", "Could not obtain the Getter object for file %s and event %s !", GetTitle(), fEventFolderName.Data()) ;   
454     return kFALSE;
455   } 
456   TString opt("Digits") ; 
457   if(gime->VersionExists(opt) ) { 
458     Error( "Init", "Give a version name different from %s", fEventFolderName.Data() ) ;
459     fInit = kFALSE ; 
460   }
461
462   fFirstEvent = 0 ; 
463   fLastEvent = fFirstEvent ; 
464   
465   if(fManager)
466     fInput = fManager->GetNinputs() ; 
467   else 
468     fInput           = 1 ;
469
470   fInputFileNames  = new TString[fInput] ; 
471   fEventNames      = new TString[fInput] ; 
472   fInputFileNames[0] = GetTitle() ; 
473   fEventNames[0]     = fEventFolderName.Data() ; 
474   Int_t index ; 
475   for (index = 1 ; index < fInput ; index++) {
476     fInputFileNames[index] = dynamic_cast<AliStream*>(fManager->GetInputStream(index))->GetFileName(0); 
477     TString tempo = fManager->GetInputFolderName(index) ;
478     fEventNames[index] = tempo.Remove(tempo.Length()-1) ; // strip of the stream number added bt fManager 
479   }
480   
481   //to prevent cleaning of this object while GetEvent is called
482   gime->EmcalLoader()->GetDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
483   
484   return fInit ;    
485 }
486
487 //____________________________________________________________________________ 
488 void AliEMCALDigitizer::InitParameters()
489 { // Tune parameters - 24-nov-04
490
491   fMeanPhotonElectron = 3300 ; // electrons per GeV 
492   fPinNoise           = 0.004; 
493   if (fPinNoise == 0. ) 
494     Warning("InitParameters", "No noise added\n") ; 
495   fDigitThreshold     = fPinNoise * 3; // 3 * sigma
496   fTimeResolution     = 0.3e-9 ; // 300 psc
497   fTimeSignalLength   = 1.0e-9 ;
498
499   fADCchannelEC    = 0.00305; // 200./65536 - width of one ADC channel in GeV
500   fADCpedestalEC   = 0.009 ;  // GeV
501   fNADCEC          = (Int_t) TMath::Power(2,16) ;  // number of channels in Tower ADC - 65536
502
503   fTimeThreshold      = 0.001*10000000 ; // Means 1 MeV in terms of SDigits amplitude ??
504   // hists. for control; no hists on default
505   fControlHists = 0;
506   fHists        = 0;
507 }
508
509 //__________________________________________________________________
510 void AliEMCALDigitizer::MixWith(TString alirunFileName, TString eventFolderName)
511 {
512   // Allows to produce digits by superimposing background and signal event.
513   // It is assumed, that headers file with SIGNAL events is opened in 
514   // the constructor. 
515   // Sets the BACKGROUND event, with which the SIGNAL event is to be mixed 
516   // Thus we avoid writing (changing) huge and expensive 
517   // backgound files: all output will be writen into SIGNAL, i.e. 
518   // opened in constructor file. 
519   //
520   // One can open as many files to mix with as one needs.
521   // However only Sdigits with the same name (i.e. constructed with the same SDigitizer)
522   // can be mixed.
523
524   if( strcmp(GetName(), "") == 0 )
525     Init() ;
526   
527   if(fManager){
528     Error("MixWith", "Cannot use this method under AliRunDigitizer") ;
529     return ;
530   } 
531   // looking for file which contains AliRun
532   if (gSystem->AccessPathName(alirunFileName)) {// file does not exist
533     Error("MixWith", "File %s does not exist!", alirunFileName.Data()) ;
534     return ; 
535   }
536   // looking for the file which contains SDigits
537   AliEMCALGetter * gime = AliEMCALGetter::Instance() ; 
538   TString fileName( gime->GetSDigitsFileName() ) ; 
539     if ( eventFolderName != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name 
540       fileName = fileName.ReplaceAll(".root", "") + "_" + eventFolderName + ".root" ;
541     if ( (gSystem->AccessPathName(fileName)) ) { 
542       Error("MixWith", "The file %s does not exist!", fileName.Data()) ;
543       return ;
544     }
545     // need to increase the arrays
546     TString tempo = fInputFileNames[fInput-1] ; 
547     delete [] fInputFileNames ; 
548     fInputFileNames = new TString[fInput+1] ; 
549     fInputFileNames[fInput-1] = tempo ; 
550  
551     tempo = fEventNames[fInput-1] ; 
552     delete [] fEventNames ; 
553     fEventNames = new TString[fInput+1] ; 
554     fEventNames[fInput-1] = tempo ; 
555
556     fInputFileNames[fInput] = alirunFileName ; 
557     fEventNames[fInput]     = eventFolderName ;
558     fInput++ ;
559 }  
560
561 void AliEMCALDigitizer::Print1(Option_t * option)
562 { // 19-nov-04 - just for convinience
563   Print(); 
564   PrintDigits(option);
565 }
566
567 //__________________________________________________________________
568 void AliEMCALDigitizer::Print()const 
569 {
570   // Print Digitizer's parameters
571   printf("Print: \n------------------- %s -------------", GetName() ) ; 
572   if( strcmp(fEventFolderName.Data(), "") != 0 ){
573     printf(" Writing Digits to branch with title  %s\n", fEventFolderName.Data()) ;
574     
575     Int_t nStreams ; 
576     if (fManager) 
577       nStreams =  GetNInputStreams() ;
578     else 
579       nStreams = fInput ; 
580     
581     Int_t index = 0 ;  
582     for (index = 0 ; index < nStreams ; index++) {  
583       TString tempo(fEventNames[index]) ; 
584       tempo += index ;
585       AliEMCALGetter * gime = AliEMCALGetter::Instance(fInputFileNames[index], tempo) ; 
586       TString fileName( gime->GetSDigitsFileName() ) ; 
587       if ( fEventNames[index] != AliConfig::GetDefaultEventFolderName()) // only if not the default folder name 
588         fileName = fileName.ReplaceAll(".root", "") + "_" + fEventNames[index]  + ".root" ;
589       printf ("Adding SDigits from %s %s\n", fInputFileNames[index].Data(), fileName.Data()) ; 
590     }
591     AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName) ; 
592     printf("\nWriting digits to %s", gime->GetDigitsFileName().Data()) ;
593     
594     printf("\nWith following parameters:\n") ;
595     
596     printf("    Electronics noise in EMC (fPinNoise) = %f\n", fPinNoise) ;
597     printf("    Threshold  in EMC  (fDigitThreshold) = %f\n", fDigitThreshold)  ;
598     printf("---------------------------------------------------\n")  ;
599   }
600   else
601     printf("Print: AliEMCALDigitizer not initialized") ; 
602 }
603
604 //__________________________________________________________________
605 void AliEMCALDigitizer::PrintDigits(Option_t * option){
606     
607   AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName) ; 
608   TClonesArray * digits  = gime->Digits() ;
609   TClonesArray * sdigits = gime->SDigits() ;
610   
611   printf("\n #Digits: %d : sdigits %d ", digits->GetEntriesFast(), sdigits->GetEntriesFast()) ; 
612   printf("\n event %d", gAlice->GetEvNumber()) ;
613   
614   if(strstr(option,"all")){  
615     //loop over digits
616     AliEMCALDigit * digit;
617     printf("\nEMC digits (with primaries):\n")  ;
618     printf("\n   Id  Amplitude    Time          Index Nprim: Primaries list \n") ;    
619     Int_t index ;
620     for (index = 0 ; index < digits->GetEntries() ; index++) {
621       digit = dynamic_cast<AliEMCALDigit *>(digits->At(index)) ;
622       printf("\n%6d  %8d    %6.5e %4d      %2d : ",
623               digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;  
624       Int_t iprimary;
625       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
626         printf("%d ",digit->GetPrimary(iprimary+1) ) ; 
627       }
628     }   
629   }
630   printf("\n");
631 }
632
633 //__________________________________________________________________
634 Float_t AliEMCALDigitizer::TimeOfNoise(void)
635 {  
636   // Calculates the time signal generated by noise
637   //PH  Info("TimeOfNoise", "Change me") ; 
638   return gRandom->Rndm() * 1.28E-5;
639 }
640
641 //__________________________________________________________________
642 void AliEMCALDigitizer::Unload() 
643 {  
644   // Unloads the SDigits and Digits
645   Int_t i ; 
646   for(i = 1 ; i < fInput ; i++){
647     TString tempo(fEventNames[i]) ; 
648     tempo += i ;
649     AliEMCALGetter * gime = AliEMCALGetter::Instance(fInputFileNames[i], tempo) ; 
650     gime->EmcalLoader()->UnloadSDigits() ; 
651   }
652   
653   AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName) ; 
654   gime->EmcalLoader()->UnloadDigits() ; 
655 }
656
657 //_________________________________________________________________________________________
658 void AliEMCALDigitizer::WriteDigits()
659 {
660
661   // Makes TreeD in the output file. 
662   // Check if branch already exists: 
663   //   if yes, exit without writing: ROOT TTree does not support overwriting/updating of 
664   //      already existing branches. 
665   //   else creates branch with Digits, named "EMCAL", title "...",
666   //      and branch "AliEMCALDigitizer", with the same title to keep all the parameters
667   //      and names of files, from which digits are made.
668
669   AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle()) ; 
670
671   const TClonesArray * digits = gime->Digits() ; 
672   TTree * treeD = gime->TreeD(); 
673
674   // -- create Digits branch
675   Int_t bufferSize = 32000 ;    
676   TBranch * digitsBranch = treeD->Branch("EMCAL","TClonesArray",&digits,bufferSize);
677   digitsBranch->SetTitle(fEventFolderName);
678   digitsBranch->Fill() ;
679   
680   gime->WriteDigits("OVERWRITE");
681   gime->WriteDigitizer("OVERWRITE");
682
683   Unload() ; 
684
685 }
686
687 void AliEMCALDigitizer::Browse(TBrowser* b)
688 {
689   if(fHists) b->Add(fHists);
690   TTask::Browse(b);
691 }
692
693 TList *AliEMCALDigitizer::BookControlHists(const int var)
694 { // 22-nov-04
695   Info("BookControlHists"," started ");
696   gROOT->cd();
697   AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName);
698   const AliEMCALGeometry *geom = gime->EMCALGeometry() ;
699   if(var>=1){
700     new TH1F("hDigiN",  "#EMCAL digits with fAmp > fDigitThreshold", 
701     fNADCEC+1, -0.5, Double_t(fNADCEC));
702     new TH1F("HDigiSumEnergy","Sum.EMCAL energy from digi", 1000, 0.0, 200.);
703     new TH1F("hDigiAmp",  "EMCAL digital amplitude", fNADCEC+1, -0.5, Double_t(fNADCEC));
704     new TH1F("hDigiEnergy","EMCAL cell energy", 2000, 0.0, 200.);
705     new TH1F("hDigiAbsId","EMCAL absId cells with fAmp > fDigitThreshold ",
706     geom->GetNCells(), 0.5, Double_t(geom->GetNCells())+0.5);
707   }
708   fHists = sv::MoveHistsToList("EmcalDigiControlHists", kFALSE);
709   return fHists;
710 }
711
712 void AliEMCALDigitizer::SaveHists(const char* name, Bool_t kSingleKey, const char* opt)
713 {
714   sv::SaveListOfHists(fHists, name, kSingleKey, opt); 
715 }