1 /*************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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.
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
44 //*-- Author : Sahal Yacoob (LBL)
45 // based on : AliPHOSSDigitzer
46 //////////////////////////////////////////////////////////////////////////////
49 // --- ROOT system ---
56 #include "TBenchmark.h"
57 // --- Standard library ---
60 // --- AliRoot header files ---
62 #include "AliEMCALDigit.h"
63 #include "AliEMCALHit.h"
64 #include "AliEMCALSDigitizer.h"
65 #include "AliEMCALGeometry.h"
66 #include "AliEMCALv1.h"
67 #include "AliEMCALGetter.h"
69 ClassImp(AliEMCALSDigitizer)
72 //____________________________________________________________________________
73 AliEMCALSDigitizer::AliEMCALSDigitizer():TTask("AliEMCALSDigitizer","")
76 fA = fB = fNevents = 0 ;
77 fTowerPrimThreshold = fPreShowerPrimThreshold = fPhotonElectronFactor = 0. ;
78 fHits = fSDigits = fSDigits = 0 ;
80 fIsInitialized = kFALSE ;
85 //____________________________________________________________________________
86 AliEMCALSDigitizer::AliEMCALSDigitizer(const char* headerFile, const char *sDigitsTitle):TTask(sDigitsTitle, headerFile)
91 fTowerPrimThreshold = 0.01 ;
92 fPreShowerPrimThreshold = 0.0001 ;
94 fPhotonElectronFactor = 5000. ; // photoelectrons per GeV
98 //____________________________________________________________________________
99 AliEMCALSDigitizer::~AliEMCALSDigitizer()
107 //____________________________________________________________________________
108 void AliEMCALSDigitizer::Init(){
110 // Initialization: open root-file, allocate arrays for hits and sdigits,
111 // attach task SDigitizer to the list of PHOS tasks
113 // Initialization can not be done in the default constructor
114 //============================================================= YS
115 // The initialisation is now done by the getter
117 if( strcmp(GetTitle(), "") == 0 )
118 SetTitle("galice.root") ;
120 AliEMCALGetter * gime = AliEMCALGetter::GetInstance(GetTitle(), GetName(), "update") ;
122 cerr << "ERROR: AliEMCALSDigitizer::Init -> Could not obtain the Getter object !" << endl ;
126 gime->PostSDigits( GetName(), GetTitle() ) ;
128 TString sdname(GetName() );
130 sdname.Append(GetTitle() ) ;
131 SetName(sdname.Data()) ;
132 gime->PostSDigitizer(this) ;
136 //____________________________________________________________________________
137 void AliEMCALSDigitizer::Exec(Option_t *option) {
139 // Collects all hits in the same active volume into digit
141 if( strcmp(GetName(), "") == 0 )
144 if (strstr(option, "print") ) {
149 if(strstr(option,"tim"))
150 gBenchmark->Start("EMCALSDigitizer");
153 //Check, if this branch already exits
154 gAlice->GetEvent(0) ;
156 TString sdname(GetName()) ;
157 sdname.Remove(sdname.Index(GetTitle())-1) ;
159 if(gAlice->TreeS() ) {
160 TObjArray * lob = static_cast<TObjArray*>(gAlice->TreeS()->GetListOfBranches()) ;
162 TBranch * branch = 0 ;
163 Bool_t emcalfound = kFALSE, sdigitizerfound = kFALSE ;
165 while ( (branch = (static_cast<TBranch*>(next()))) && (!emcalfound || !sdigitizerfound) ) {
166 if ( (strcmp(branch->GetName(), "EMCAL")==0) && (strcmp(branch->GetTitle(), GetName())==0) )
169 else if ( (strcmp(branch->GetName(), "AliEMCALSDigitizer")==0) && (strcmp(branch->GetTitle(), sdname)==0) )
170 sdigitizerfound = kTRUE ;
173 if ( emcalfound || sdigitizerfound ) {
174 cerr << "WARNING: AliEMCALSDigitizer::Exec -> SDigits and/or SDigitizer branch with name " << GetName()
175 << " already exits" << endl ;
180 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
181 Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
184 for(ievent = 0; ievent < nevents; ievent++){
185 gime->Event(ievent,"H") ;
186 const TClonesArray * fHits = gime->Hits() ;
187 TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
191 //Collects all hits in the same active volume into digit
196 //Now made SDigits from hits, for EMCAL it is the same, so just copy
198 // for (itrack=0; itrack < gAlice->GetNtrack(); itrack++){
199 //gime->Track(itrack);
200 //=========== Get the EMCAL branch from Hits Tree for the Primary track itrack
203 for ( i = 0 ; i < fHits->GetEntries() ; i++ ) { // loop over all hits (hit = deposited energy/layer/entering particle)
204 AliEMCALHit * hit = dynamic_cast<AliEMCALHit*>(fHits->At(i)) ;
205 AliEMCALDigit * curSDigit = 0 ;
206 AliEMCALDigit * sdigit = 0 ;
207 Bool_t newsdigit = kTRUE;
211 // Assign primary number only if deposited energy is significant
213 if( (!hit->IsInPreShower() && hit->GetEnergy() > fTowerPrimThreshold) ||
214 (hit->IsInPreShower() && hit->GetEnergy() > fPreShowerPrimThreshold)) {
215 curSDigit = new AliEMCALDigit( hit->GetPrimary(),
216 hit->GetIparent(),Layer2TowerID(hit->GetId(),hit->IsInPreShower()),
217 Digitize(hit->GetEnergy()),
220 curSDigit = new AliEMCALDigit( -1 ,
222 Layer2TowerID(hit->GetId(),hit->IsInPreShower()),
223 Digitize(hit->GetEnergy()),
227 for(check= 0; check < nSdigits ; check++) {
228 sdigit = (AliEMCALDigit *)sdigits->At(check);
229 if( sdigit->GetId() == curSDigit->GetId()) { // Are we in the same tower or the same preshower ?
230 *sdigit = *sdigit + *curSDigit;
236 new((*sdigits)[nSdigits]) AliEMCALDigit(*curSDigit);
239 } // loop over all hits (hit = deposited energy/layer/entering particle)
240 // } // loop over tracks
244 nSdigits = sdigits->GetEntriesFast() ;
245 sdigits->Expand(nSdigits) ;
248 const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
250 Int_t lastPreShowerIndex = nSdigits - 1 ;
251 if (!(dynamic_cast<AliEMCALDigit *>(sdigits->At(lastPreShowerIndex))->IsInPreShower()))
252 lastPreShowerIndex = -2;
253 Int_t firstPreShowerIndex = 100000 ;
255 AliEMCALDigit * sdigit = 0 ;
256 for ( index = 0; index < nSdigits ; index++) {
257 sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(index) ) ;
258 if (sdigit->IsInPreShower() ){
259 firstPreShowerIndex = index ;
264 AliEMCALDigit * preshower ;
265 AliEMCALDigit * tower ;
266 Int_t lastIndex = lastPreShowerIndex +1 ;
269 for (index = firstPreShowerIndex ; index <= lastPreShowerIndex; index++) {
270 preshower = dynamic_cast<AliEMCALDigit *>(sdigits->At(index) );
271 Bool_t towerFound = kFALSE ;
273 for (jndex = 0; jndex < firstPreShowerIndex; jndex++) {
274 tower = dynamic_cast<AliEMCALDigit *>(sdigits->At(jndex) );
275 if ( (preshower->GetId() - (geom->GetNZ() * geom->GetNPhi()) ) == tower->GetId() ) {
276 Float_t towerEnergy = static_cast<Float_t>(tower->GetAmp()) ;
277 Float_t preshoEnergy = static_cast<Float_t>(preshower->GetAmp()) ;
278 towerEnergy +=preshoEnergy ;
279 *tower = *tower + *preshower ; // and add preshower multiplied by layer ratio to tower
280 tower->SetAmp(static_cast<Int_t>(TMath::Ceil(towerEnergy))) ;
286 new((*sdigits)[lastIndex]) AliEMCALDigit(*preshower);
287 AliEMCALDigit * temp = dynamic_cast<AliEMCALDigit *>(sdigits->At(lastIndex)) ;
288 temp->SetId(temp->GetId() - (geom->GetNZ() * geom->GetNPhi()) ) ;
294 Int_t NPrimarymax = -1 ;
295 for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) {
296 sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(i)) ;
297 sdigit->SetIndexInList(i) ;
300 for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) {
301 if (((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) > NPrimarymax)
302 NPrimarymax = ((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) ;
304 if(gAlice->TreeS() == 0)
305 gAlice->MakeTree("S") ;
307 //Make (if necessary) branches
309 if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
310 file = new char[strlen(gAlice->GetBaseFile())+20] ;
311 sprintf(file,"%s/EMCAL.SDigits.root",gAlice->GetBaseFile()) ;
314 TDirectory *cwd = gDirectory;
316 //First list of sdigits
317 Int_t bufferSize = 32000 ;
318 TBranch * sdigitsBranch = gAlice->TreeS()->Branch("EMCAL",&sdigits,bufferSize);
319 sdigitsBranch->SetTitle(sdname);
322 sdigitsBranch->SetFile(file);
323 TIter next( sdigitsBranch->GetListOfBranches());
325 while ((subbr=(TBranch*)next())) {
326 subbr->SetFile(file);
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);
338 sdigitizerBranch->SetFile(file);
339 TIter next( sdigitizerBranch->GetListOfBranches());
341 while ((subbr=(TBranch*)next())) {
342 subbr->SetFile(file);
348 gAlice->TreeS()->Fill() ;
349 gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
351 if(strstr(option,"deb"))
352 PrintSDigits(option) ;
356 if(strstr(option,"tim")){
357 gBenchmark->Stop("EMCALSDigitizer");
358 cout << "AliEMCALSDigitizer:" << endl ;
359 cout << " took " << gBenchmark->GetCpuTime("EMCALSDigitizer") << " seconds for SDigitizing "
360 << gBenchmark->GetCpuTime("EMCALSDigitizer")/fNevents << " seconds per event " << endl ;
366 //__________________________________________________________________
367 void AliEMCALSDigitizer::SetSDigitsBranch(const char * title ){
369 // Setting title to branch SDigits
371 TString stitle(title) ;
373 // check if branch with title already exists
374 TBranch * sdigitsBranch =
375 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("EMCAL")) ;
376 TBranch * sdigitizerBranch =
377 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("AliEMCALSDigitizer")) ;
378 const char * sdigitsTitle = sdigitsBranch ->GetTitle() ;
379 const char * sdigitizerTitle = sdigitizerBranch ->GetTitle() ;
380 if ( stitle.CompareTo(sdigitsTitle)==0 || stitle.CompareTo(sdigitizerTitle)==0 ){
381 cerr << "ERROR: AliEMCALSdigitizer::SetSDigitsBranch -> Cannot overwrite existing branch with title " << title << endl ;
385 cout << "AliEMCALSdigitizer::SetSDigitsBranch -> Changing SDigits file from " << GetName() << " to " << title << endl ;
389 // Post to the WhiteBoard
390 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
391 gime->PostSDigits( title, GetTitle()) ;
395 //__________________________________________________________________
396 void AliEMCALSDigitizer::Print(Option_t* option)const{
397 cout << "------------------- "<< GetName() << " -------------" << endl ;
398 cout << " Writing SDigitis to branch with title " << GetName() << endl ;
399 cout << " with digitization parameters A = " << fA << endl ;
400 cout << " B = " << fB << endl ;
401 cout << " Threshold for Primary assignment in Tower = " << fTowerPrimThreshold << endl ;
402 cout << " Threshold for Primary assignment in PreShower = " << fPreShowerPrimThreshold << endl ;
403 cout << "---------------------------------------------------"<<endl ;
406 //__________________________________________________________________
407 Bool_t AliEMCALSDigitizer::operator==( AliEMCALSDigitizer const &sd )const{
408 if( (fA==sd.fA)&&(fB==sd.fB)&&
409 (fTowerPrimThreshold==sd.fTowerPrimThreshold) &&
410 (fPreShowerPrimThreshold==sd.fPreShowerPrimThreshold))
415 //__________________________________________________________________
416 void AliEMCALSDigitizer::PrintSDigits(Option_t * option){
417 //Prints list of digits produced at the current pass of AliEMCALDigitizer
419 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
420 TString sdname(GetName()) ;
421 sdname.Remove(sdname.Index(GetTitle())-1) ;
422 const TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
424 cout << "AliEMCALSDigitiser: event " << gAlice->GetEvNumber() << endl ;
425 cout << " Number of entries in SDigits list " << sdigits->GetEntriesFast() << endl ;
427 if(strstr(option,"all")||strstr(option,"EMC")){
430 AliEMCALDigit * digit;
431 cout << "SDigit Id " << " Amplitude " << " Time " << " Index " << " Nprim " << " Primaries list " << endl;
433 for (index = 0 ; index < sdigits->GetEntries() ; index++) {
434 digit = (AliEMCALDigit * ) sdigits->At(index) ;
435 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " " << digit->GetTime()
436 << setw(6) << digit->GetIndexInList() << " "
437 << setw(5) << digit->GetNprimary() <<" ";
440 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
441 cout << " " << digit->GetPrimary(iprimary+1) << " ";
447 //________________________________________________________________________
448 Int_t AliEMCALSDigitizer::Layer2TowerID(Int_t ihit, Bool_t preshower){
449 // Method to Transform from Hit Id to Digit Id
450 // This function should be one to one
451 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
452 const AliEMCALGeometry * geom = gime->EMCALGeometry();
453 Int_t ieta = ((ihit-1)/geom->GetNPhi())%geom->GetNZ(); // eta Tower Index
454 Int_t iphi = (ihit-1)%(geom->GetNPhi())+1; //phi Tower Index
458 if (preshower)ipre = 1;
459 if (iphi > 0 && ieta >= 0){
460 it = iphi+ieta*geom->GetNPhi() + ipre*geom->GetNPhi()*geom->GetNZ();
463 cerr << " AliEMCALSDigitizer::Layer2TowerID() -- there is an error "<< endl << "Eta number = "
464 << ieta << "Phi number = " << iphi << endl ;
466 } // end if iphi>0 && ieta>0
468 //_______________________________________________________________________________________
469 // void AliEMCALSDigitizer::TestTowerID(void)
473 // Bool_t preshower = kFALSE;
474 // for (j = 0 ; j < 10 ; j++){ // loop over hit id
476 // for (i = 0 ; i <= 2 ; i++){ // loop over
477 // Int_t k = i*96*144+j*144+1;
478 // cout << " Hit Index = " << k << " " << j*10 << " TOWERID = " << Layer2TowerID(k, preshower) << endl ;