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