]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALSDigitizer.cxx
SDigits has become a copy of hits and threshold to associate primary particles has...
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALSDigitizer.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 // This is a TTask that makes SDigits out of Hits
20 // A Summable Digits is the sum of all hits originating 
21 // from one in one tower of the EMCAL 
22 // A threshold for assignment of the primary to SDigit is applied 
23 // SDigits are written to TreeS, branch "EMCAL"
24 // AliEMCALSDigitizer with all current parameters is written 
25 // to TreeS branch "AliEMCALSDigitizer".
26 // Both branches have the same title. If necessary one can produce 
27 // another set of SDigits with different parameters. Two versions
28 // can be distunguished using titles of the branches.
29 // User case:
30 // root [0] AliEMCALSDigitizer * s = new AliEMCALSDigitizer("galice.root")
31 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
32 // root [1] s->ExecuteTask()
33 //             // Makes SDigitis for all events stored in galice.root
34 // root [2] s->SetPedestalParameter(0.001)
35 //             // One can change parameters of digitization
36 // root [3] s->SetSDigitsBranch("Redestal 0.001")
37 //             // and write them into the new branch
38 // root [4] s->ExeciteTask("deb all tim")
39 //             // available parameters:
40 //             deb - print # of produced SDigitis
41 //             deb all  - print # and list of produced SDigits
42 //             tim - print benchmarking information
43 //
44 //*-- Author : Sahal Yacoob (LBL)
45 // based on  : AliPHOSSDigitzer 
46 //////////////////////////////////////////////////////////////////////////////
47
48
49 // --- ROOT system ---
50 #include "TFile.h"
51 #include "TTask.h"
52 #include "TTree.h"
53 #include "TSystem.h"
54 #include "TROOT.h"
55 #include "TFolder.h"
56 #include "TBenchmark.h"
57 // --- Standard library ---
58 #include <iomanip.h>
59
60 // --- AliRoot header files ---
61 #include "AliRun.h"
62 #include "AliEMCALDigit.h"
63 #include "AliEMCALHit.h"
64 #include "AliEMCALSDigitizer.h"
65 #include "AliEMCALGeometry.h"
66 #include "AliEMCALv1.h"
67 #include "AliEMCALGetter.h"
68
69 ClassImp(AliEMCALSDigitizer)
70
71            
72 //____________________________________________________________________________ 
73   AliEMCALSDigitizer::AliEMCALSDigitizer():TTask("AliEMCALSDigitizer","") 
74 {
75   // ctor
76   fA = fB =  fNevents = 0 ; 
77   fTowerPrimThreshold = fPreShowerPrimThreshold = fPhotonElectronFactor  = 0. ;
78   fHits = fSDigits = fSDigits = 0 ;
79
80   fIsInitialized = kFALSE ;
81
82
83 }
84
85 //____________________________________________________________________________ 
86 AliEMCALSDigitizer::AliEMCALSDigitizer(const char* headerFile, const char *sDigitsTitle):TTask(sDigitsTitle, headerFile)
87 {
88   // ctor
89   fA = 0;
90   fB = 10000000.;
91   fTowerPrimThreshold = 0.01 ;
92   fPreShowerPrimThreshold = 0.0001 ; 
93   fNevents = 0 ; 
94   fPhotonElectronFactor = 5000. ; // photoelectrons per GeV 
95   Init();
96 }
97
98 //____________________________________________________________________________ 
99 AliEMCALSDigitizer::~AliEMCALSDigitizer()
100 {
101   // dtor
102   if(fSDigits)
103     delete fSDigits ;
104   if(fHits)
105     delete fHits ;
106 }
107 //____________________________________________________________________________ 
108 void AliEMCALSDigitizer::Init(){
109
110   // Initialization: open root-file, allocate arrays for hits and sdigits,
111   // attach task SDigitizer to the list of PHOS tasks
112   // 
113   // Initialization can not be done in the default constructor
114   //============================================================= YS
115   //  The initialisation is now done by the getter
116
117 if( strcmp(GetTitle(), "") == 0 )
118     SetTitle("galice.root") ;
119   
120   AliEMCALGetter * gime = AliEMCALGetter::GetInstance(GetTitle(), GetName(), "update") ;     
121   if ( gime == 0 ) {
122     cerr << "ERROR: AliEMCALSDigitizer::Init -> Could not obtain the Getter object !" << endl ; 
123     return ;
124   } 
125   
126   gime->PostSDigits( GetName(), GetTitle() ) ; 
127
128   TString sdname(GetName() );
129   sdname.Append(":") ;
130   sdname.Append(GetTitle() ) ;
131   SetName(sdname) ;
132   gime->PostSDigitizer(this) ;
133  
134   // create a folder on the white board //YSAlice/WhiteBoard/SDigits/EMCAL/headerFile/sdigitsTitle
135  
136  }
137 //____________________________________________________________________________
138 void AliEMCALSDigitizer::Exec(Option_t *option) { 
139
140 // Collects all hits in the same active volume into digit
141  
142   if( strcmp(GetName(), "") == 0 )
143     Init() ;
144   
145   if (strstr(option, "print") ) {
146     Print("") ; 
147     return ; 
148   }
149   
150   if(strstr(option,"tim"))
151     gBenchmark->Start("EMCALSDigitizer");
152   //Check, if this branch already exits
153   gAlice->GetEvent(0) ;
154   if(gAlice->TreeS() ) {
155     TObjArray * lob = static_cast<TObjArray*>(gAlice->TreeS()->GetListOfBranches()) ;
156     TIter next(lob) ; 
157     TBranch * branch = 0 ;  
158     Bool_t emcalfound = kFALSE, sdigitizerfound = kFALSE ; 
159     
160     while ( (branch = (static_cast<TBranch*>(next()))) && (!emcalfound || !sdigitizerfound) ) {
161       if ( (strcmp(branch->GetName(), "EMCAL")==0) && (strcmp(branch->GetTitle(), GetName())==0) ) 
162         emcalfound = kTRUE ;
163       
164       else if ( (strcmp(branch->GetName(), "AliEMCALSDigitizer")==0) && (strcmp(branch->GetTitle(), GetName())==0) ) 
165         sdigitizerfound = kTRUE ; 
166     }
167     
168     if ( emcalfound || sdigitizerfound ) {
169       cerr << "WARNING: AliEMCALSDigitizer::Exec -> SDigits and/or SDigitizer branch with name " << GetName() 
170            << " already exits" << endl ;
171       return ; 
172     }   
173   }  
174   TString sdname(GetName()) ;
175   sdname.Remove(sdname.Index(GetTitle())-1) ;
176   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;  
177   Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ; 
178   
179   Int_t ievent ;
180   for(ievent = 0; ievent < nevents; ievent++){
181     gime->Event(ievent,"H") ;
182     const TClonesArray * fHits = gime->Hits() ;
183     TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
184     sdigits->Clear();
185     Int_t nSdigits = 0 ;
186
187   //Collects all hits in the same active volume into digit
188   
189   
190     
191     
192     //Now made SDigits from hits, for EMCAL it is the same, so just copy    
193     //  Int_t itrack ;
194     // for (itrack=0; itrack < gAlice->GetNtrack(); itrack++){
195     //gime->Track(itrack);  
196       //=========== Get the EMCAL branch from Hits Tree for the Primary track itrack
197      
198       Int_t i;
199       for ( i = 0 ; i < fHits->GetEntries() ; i++ ) { // loop over all hits (hit = deposited energy/layer/entering particle)
200         AliEMCALHit * hit = dynamic_cast<AliEMCALHit*>(fHits->At(i)) ;
201         AliEMCALDigit * curSDigit = 0 ;
202         AliEMCALDigit * sdigit = 0 ;
203         Bool_t newsdigit = kTRUE; 
204
205         
206
207         // Assign primary number only if deposited energy is significant
208           
209         if( (!hit->IsInPreShower() && hit->GetEnergy() > fTowerPrimThreshold) || 
210             (hit->IsInPreShower() && hit->GetEnergy() > fPreShowerPrimThreshold)) {
211           curSDigit =  new AliEMCALDigit( hit->GetPrimary(),
212                                           hit->GetIparent(),Layer2TowerID(hit->GetId(),hit->IsInPreShower()), 
213                                           Digitize(hit->GetEnergy()), 
214                                           hit->GetTime()) ;
215         } else {
216           curSDigit =  new AliEMCALDigit( -1               , 
217                                           -1               ,
218                                           Layer2TowerID(hit->GetId(),hit->IsInPreShower()), 
219                                           Digitize(hit->GetEnergy()), 
220                                           hit->GetTime() ) ;    
221         }
222         Int_t check = 0 ;
223         for(check= 0; check < nSdigits ; check++) {
224           sdigit = (AliEMCALDigit *)sdigits->At(check);
225           if( sdigit->GetId() == curSDigit->GetId()) { // Are we in the same tower or the same preshower ?              
226               *sdigit = *sdigit + *curSDigit;
227    
228               newsdigit = kFALSE;
229           }
230         }
231         if (newsdigit) { 
232           new((*sdigits)[nSdigits])  AliEMCALDigit(*curSDigit);
233           nSdigits++ ;  
234         }
235       }  // loop over all hits (hit = deposited energy/layer/entering particle)
236       //    } // loop over tracks
237       
238       sdigits->Sort() ;
239       
240       nSdigits = sdigits->GetEntriesFast() ;
241       sdigits->Expand(nSdigits) ;
242       
243       //  Int_t i ;
244      const AliEMCALGeometry * geom = gime->EMCALGeometry() ; 
245      
246      Int_t lastPreShowerIndex = nSdigits - 1 ;
247      if (!(dynamic_cast<AliEMCALDigit *>(sdigits->At(lastPreShowerIndex))->IsInPreShower()))
248        lastPreShowerIndex = -2; 
249      Int_t firstPreShowerIndex = 100000 ; 
250      Int_t index ; 
251      AliEMCALDigit * sdigit = 0 ;
252      for ( index = 0; index < nSdigits ; index++) {
253        sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(index) ) ;
254        if (sdigit->IsInPreShower() ){ 
255          firstPreShowerIndex = index ;
256          break ;
257        }
258      }
259      
260     AliEMCALDigit * preshower ;
261     AliEMCALDigit * tower ;
262     Int_t lastIndex = lastPreShowerIndex +1 ; 
263
264
265     for (index = firstPreShowerIndex ; index <= lastPreShowerIndex; index++) {
266       preshower = dynamic_cast<AliEMCALDigit *>(sdigits->At(index) ); 
267       Bool_t towerFound = kFALSE ;
268       Int_t jndex ; 
269       for (jndex = 0; jndex < firstPreShowerIndex; jndex++) {
270         tower  = dynamic_cast<AliEMCALDigit *>(sdigits->At(jndex) ); 
271         if ( (preshower->GetId() - (geom->GetNZ() * geom->GetNPhi()) ) == tower->GetId() ) {      
272           Float_t towerEnergy  = static_cast<Float_t>(tower->GetAmp()) ; 
273           Float_t preshoEnergy = static_cast<Float_t>(preshower->GetAmp()) ; 
274           towerEnergy +=preshoEnergy ; 
275           *tower = *tower + *preshower    ; // and add preshower multiplied by layer ratio to tower
276           tower->SetAmp(static_cast<Int_t>(TMath::Ceil(towerEnergy))) ; 
277           towerFound = kTRUE ;
278         }
279       }
280       if ( !towerFound ) { 
281
282         new((*sdigits)[lastIndex])  AliEMCALDigit(*preshower);
283         AliEMCALDigit * temp = dynamic_cast<AliEMCALDigit *>(sdigits->At(lastIndex)) ;
284         temp->SetId(temp->GetId() - (geom->GetNZ() * geom->GetNPhi()) ) ;
285         lastIndex++ ; 
286       }
287     }
288     
289     sdigits->Sort() ;
290     Int_t NPrimarymax = -1 ; 
291     for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) { 
292       sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(i)) ;
293       sdigit->SetIndexInList(i) ;
294     }
295     
296     for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) { 
297          if (((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) > NPrimarymax)
298            NPrimarymax = ((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) ;
299        }
300        if(gAlice->TreeS() == 0)
301         gAlice->MakeTree("S") ;
302       cout << " AliEMCALSDigitizer:: NPrimaryMax = " << NPrimarymax << endl;
303       
304       
305       
306       //Make (if necessary) branches    
307       char * file =0;
308       if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
309         file = new char[strlen(gAlice->GetBaseFile())+20] ;
310         sprintf(file,"%s/EMCAL.SDigits.root",gAlice->GetBaseFile()) ;
311       }
312       
313       TDirectory *cwd = gDirectory;
314       
315       //First list of sdigits
316       Int_t bufferSize = 32000 ;    
317       TBranch * sdigitsBranch = gAlice->TreeS()->Branch("EMCAL",&sdigits,bufferSize);
318       sdigitsBranch->SetTitle(sdname);
319       cout << " AliEMCALSDigitizer::Exec sdname " << sdname << endl ; 
320       
321       if (file) {
322         sdigitsBranch->SetFile(file);
323         TIter next( sdigitsBranch->GetListOfBranches());
324         TBranch * subbr;
325         while ((subbr=(TBranch*)next())) {
326           subbr->SetFile(file);
327         }   
328         cwd->cd();
329       } 
330       
331       //second - SDigitizer
332       Int_t splitlevel = 0 ;
333       AliEMCALSDigitizer * sd = this ;
334       TBranch * sdigitizerBranch = gAlice->TreeS()->Branch("AliEMCALSDigitizer","AliEMCALSDigitizer",
335                                                            &sd,bufferSize,splitlevel); 
336       sdigitizerBranch->SetTitle(sdname);
337       if (file) {
338         sdigitizerBranch->SetFile(file);
339         TIter next( sdigitizerBranch->GetListOfBranches());
340         TBranch * subbr ;
341         while ((subbr=(TBranch*)next())) {
342           subbr->SetFile(file);
343         }   
344         cwd->cd();
345         delete file;
346       }
347       
348       sdigitsBranch->Fill() ;
349       sdigitizerBranch->Fill() ;
350       gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
351       
352       if(strstr(option,"deb"))
353         PrintSDigits(option) ;
354       
355   }
356   
357   if(strstr(option,"tim")){
358     gBenchmark->Stop("EMCALSDigitizer");
359     cout << "AliEMCALSDigitizer:" << endl ;
360     cout << "   took " << gBenchmark->GetCpuTime("EMCALSDigitizer") << " seconds for SDigitizing " 
361          <<  gBenchmark->GetCpuTime("EMCALSDigitizer")/fNevents << " seconds per event " << endl ;
362     cout << endl ;
363   }
364   
365   
366 }
367 //__________________________________________________________________
368 void AliEMCALSDigitizer::SetSDigitsBranch(const char * title ){
369  
370   // Setting title to branch SDigits 
371
372   TString stitle(title) ;
373
374   // check if branch with title already exists
375   TBranch * sdigitsBranch    = 
376     static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("EMCAL")) ; 
377   TBranch * sdigitizerBranch =  
378     static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("AliEMCALSDigitizer")) ;
379   const char * sdigitsTitle    = sdigitsBranch ->GetTitle() ;  
380   const char * sdigitizerTitle = sdigitizerBranch ->GetTitle() ;
381   if ( stitle.CompareTo(sdigitsTitle)==0 || stitle.CompareTo(sdigitizerTitle)==0 ){
382     cerr << "ERROR: AliEMCALSdigitizer::SetSDigitsBranch -> Cannot overwrite existing branch with title " << title << endl ;
383     return ;
384   }
385   
386   cout << "AliEMCALSdigitizer::SetSDigitsBranch -> Changing SDigits file from " << GetName() << " to " << title << endl ;
387
388   SetName(title) ; 
389     
390   // Post to the WhiteBoard
391   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
392   gime->PostSDigits( title, GetTitle()) ; 
393
394
395 }
396 //__________________________________________________________________
397 void AliEMCALSDigitizer::Print(Option_t* option)const{
398   cout << "------------------- "<< GetName() << " -------------" << endl ;
399   cout << "   Writing SDigitis to branch with title  " << GetName() << endl ;
400   cout << "   with digitization parameters  A               = " << fA << endl ;
401   cout << "                                 B               = " << fB << endl ;
402   cout << "   Threshold for Primary assignment in Tower     = " << fTowerPrimThreshold << endl ; 
403   cout << "   Threshold for Primary assignment in PreShower = " << fPreShowerPrimThreshold << endl ; 
404   cout << "---------------------------------------------------"<<endl ;
405   
406 }
407 //__________________________________________________________________
408 Bool_t AliEMCALSDigitizer::operator==( AliEMCALSDigitizer const &sd )const{
409   if( (fA==sd.fA)&&(fB==sd.fB)&&
410       (fTowerPrimThreshold==sd.fTowerPrimThreshold) &&
411       (fPreShowerPrimThreshold==sd.fPreShowerPrimThreshold))
412     return kTRUE ;
413   else
414     return kFALSE ;
415 }
416 //__________________________________________________________________
417 void AliEMCALSDigitizer::PrintSDigits(Option_t * option){
418   //Prints list of digits produced at the current pass of AliEMCALDigitizer
419   
420   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
421   TString sdname(GetName()) ;
422   sdname.Remove(sdname.Index(GetTitle())-1) ;
423   const TClonesArray * sdigits = gime->SDigits(sdname.Data()) ; 
424   
425   cout << "AliEMCALSDigitiser: event " << gAlice->GetEvNumber() << endl ;
426   cout << "      Number of entries in SDigits list " << sdigits->GetEntriesFast() << endl ;
427   cout << endl ;
428   if(strstr(option,"all")||strstr(option,"EMC")){
429
430     //loop over digits
431     AliEMCALDigit * digit;
432     cout << "SDigit Id " << " Amplitude " <<  "     Time " << "     Index "  <<  " Nprim " << " Primaries list " <<  endl;    
433     Int_t index ;
434     for (index = 0 ; index < sdigits->GetEntries() ; index++) {
435       digit = (AliEMCALDigit * )  sdigits->At(index) ;
436       cout << setw(6)  <<  digit->GetId() << "   "  <<  setw(10)  <<  digit->GetAmp() <<   "    "  << digit->GetTime()
437            << setw(6)  <<  digit->GetIndexInList() << "    "   
438            << setw(5)  <<  digit->GetNprimary() <<"    ";
439       
440       Int_t iprimary;
441       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
442         cout << " "  <<  digit->GetPrimary(iprimary+1) << "  ";
443       cout << endl;      
444     }
445     cout <<endl;
446   }
447 }
448 //________________________________________________________________________
449 Int_t AliEMCALSDigitizer::Layer2TowerID(Int_t ihit, Bool_t preshower){
450   // Method to Transform from Hit Id to Digit Id
451   // This function should be one to one
452   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
453   const AliEMCALGeometry * geom = gime->EMCALGeometry();
454   Int_t ieta  = (ihit/geom->GetNPhi())%geom->GetNZ(); // eta Tower Index
455   Int_t iphi = (ihit-1)%(geom->GetNPhi())+1; //phi Tower Index
456   Int_t it = -10;
457   Int_t ipre = 0;
458
459   if (preshower)ipre = 1;
460   if (iphi > 0 && ieta >= 0){
461     it = iphi+ieta*geom->GetNPhi() + ipre*geom->GetNPhi()*geom->GetNZ();
462     return it;
463   }else{
464     cerr << " AliEMCALSDigitizer::Layer2TowerID() -- there is an error "<< endl << "Eta number = "
465          << ieta << "Phi number = " << iphi << endl ;
466     return it;
467   } // end if iphi>0 && ieta>0
468 }
469 //_______________________________________________________________________________________
470 // void AliEMCALSDigitizer::TestTowerID(void)
471 // {
472 //   Int_t j;
473
474 //   Bool_t preshower = kFALSE;
475 //   for (j = 0 ; j < 10 ; j++){  // loop over hit id
476 //     Int_t i;
477 //    for (i = 0 ; i <= 2 ; i++){  // loop over 
478 //      Int_t k = i*96*144+j*144+1;
479 //       cout << " Hit Index = " << k << "   " << j*10 << "   TOWERID = " <<  Layer2TowerID(k, preshower) << endl ;
480 //     }
481 //   }
482 // }