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