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 "TGeometry.h"
57 #include "TBenchmark.h"
58 // --- Standard library ---
61 // --- AliRoot header files ---
63 #include "AliHeader.h"
64 #include "AliEMCALDigit.h"
65 #include "AliEMCALHit.h"
66 #include "AliEMCALSDigitizer.h"
67 #include "AliEMCALGeometry.h"
68 #include "AliEMCALv1.h"
69 #include "AliEMCALGetter.h"
71 ClassImp(AliEMCALSDigitizer)
74 //____________________________________________________________________________
75 AliEMCALSDigitizer::AliEMCALSDigitizer():TTask("AliEMCALSDigitizer","")
78 fA = fB = fNevents = 0 ;
79 fTowerPrimThreshold = fPreShowerPrimThreshold = fPhotonElectronFactor = 0. ;
80 fHits = fSDigits = fSDigits = 0 ;
82 fIsInitialized = kFALSE ;
87 //____________________________________________________________________________
88 AliEMCALSDigitizer::AliEMCALSDigitizer(const char* headerFile, const char *sDigitsTitle):TTask(sDigitsTitle, headerFile)
93 fTowerPrimThreshold = 0.01 ;
94 fPreShowerPrimThreshold = 0.0001 ;
96 fPhotonElectronFactor = 5000. ; // photoelectrons per GeV
101 //____________________________________________________________________________
102 AliEMCALSDigitizer::~AliEMCALSDigitizer()
107 if ( fSplitFile->IsOpen() )
108 fSplitFile->Close() ;
112 //____________________________________________________________________________
113 void AliEMCALSDigitizer::Init(){
115 // Initialization: open root-file, allocate arrays for hits and sdigits,
116 // attach task SDigitizer to the list of EMCAL tasks
118 // Initialization can not be done in the default constructor
119 //============================================================= YS
120 // The initialisation is now done by the getter
122 if( strcmp(GetTitle(), "") == 0 )
123 SetTitle("galice.root") ;
125 AliEMCALGetter * gime = AliEMCALGetter::GetInstance(GetTitle(), GetName(), "update") ;
127 cerr << "ERROR: AliEMCALSDigitizer::Init -> Could not obtain the Getter object !" << endl ;
131 gime->PostSDigits( GetName(), GetTitle() ) ;
133 TString sdname(GetName() );
135 sdname.Append(GetTitle() ) ;
136 SetName(sdname.Data()) ;
137 gime->PostSDigitizer(this) ;
141 //____________________________________________________________________________
142 void AliEMCALSDigitizer::Exec(Option_t *option) {
144 // Collects all hits in the same active volume into digit
146 if( strcmp(GetName(), "") == 0 )
149 if (strstr(option, "print") ) {
154 if(strstr(option,"tim"))
155 gBenchmark->Start("EMCALSDigitizer");
158 //Check, if this branch already exits
159 gAlice->GetEvent(0) ;
161 TString sdname(GetName()) ;
162 sdname.Remove(sdname.Index(GetTitle())-1) ;
164 if(gAlice->TreeS() ) {
165 TObjArray * lob = static_cast<TObjArray*>(gAlice->TreeS()->GetListOfBranches()) ;
167 TBranch * branch = 0 ;
168 Bool_t emcalfound = kFALSE, sdigitizerfound = kFALSE ;
170 while ( (branch = (static_cast<TBranch*>(next()))) && (!emcalfound || !sdigitizerfound) ) {
171 TString thisName( GetName() ) ;
172 TString branchName( branch->GetTitle() ) ;
173 branchName.Append(":") ;
174 if ( (strcmp(branch->GetName(), "EMCAL")==0) && thisName.BeginsWith(branchName) )
177 else if ( (strcmp(branch->GetName(), "AliEMCALSDigitizer")==0) && thisName.BeginsWith(branchName) )
178 sdigitizerfound = kTRUE ;
181 if ( emcalfound || sdigitizerfound ) {
182 cerr << "WARNING: AliEMCALSDigitizer::Exec -> SDigits and/or SDigitizer branch with name " << GetName()
183 << " already exits" << endl ;
188 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
189 Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
192 for(ievent = 0; ievent < nevents; ievent++){
193 gime->Event(ievent,"H") ;
194 const TClonesArray * fHits = gime->Hits() ;
195 TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
199 //Collects all hits in the same active volume into digit
204 //Now made SDigits from hits, for EMCAL it is the same, so just copy
206 // for (itrack=0; itrack < gAlice->GetNtrack(); itrack++){
207 //gime->Track(itrack);
208 //=========== Get the EMCAL branch from Hits Tree for the Primary track itrack
211 for ( i = 0 ; i < fHits->GetEntries() ; i++ ) { // loop over all hits (hit = deposited energy/layer/entering particle)
212 AliEMCALHit * hit = dynamic_cast<AliEMCALHit*>(fHits->At(i)) ;
213 AliEMCALDigit * curSDigit = 0 ;
214 AliEMCALDigit * sdigit = 0 ;
215 Bool_t newsdigit = kTRUE;
219 // Assign primary number only if deposited energy is significant
221 if( (!hit->IsInPreShower() && hit->GetEnergy() > fTowerPrimThreshold) ||
222 (hit->IsInPreShower() && hit->GetEnergy() > fPreShowerPrimThreshold)) {
223 curSDigit = new AliEMCALDigit( hit->GetPrimary(),
224 hit->GetIparent(),Layer2TowerID(hit->GetId(),hit->IsInPreShower()),
225 Digitize(hit->GetEnergy()),
228 curSDigit = new AliEMCALDigit( -1 ,
230 Layer2TowerID(hit->GetId(),hit->IsInPreShower()),
231 Digitize(hit->GetEnergy()),
235 for(check= 0; check < nSdigits ; check++) {
236 sdigit = (AliEMCALDigit *)sdigits->At(check);
237 if( sdigit->GetId() == curSDigit->GetId()) { // Are we in the same tower or the same preshower ?
238 *sdigit = *sdigit + *curSDigit;
244 new((*sdigits)[nSdigits]) AliEMCALDigit(*curSDigit);
247 } // loop over all hits (hit = deposited energy/layer/entering particle)
248 // } // loop over tracks
252 nSdigits = sdigits->GetEntriesFast() ;
254 sdigits->Expand(nSdigits) ;
257 const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
259 Int_t lastPreShowerIndex = nSdigits - 1 ;
260 if (!(dynamic_cast<AliEMCALDigit *>(sdigits->At(lastPreShowerIndex))->IsInPreShower()))
261 lastPreShowerIndex = -2;
262 Int_t firstPreShowerIndex = 100000 ;
264 AliEMCALDigit * sdigit = 0 ;
265 for ( index = 0; index < nSdigits ; index++) {
266 sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(index) ) ;
267 if (sdigit->IsInPreShower() ){
268 firstPreShowerIndex = index ;
273 AliEMCALDigit * preshower ;
274 AliEMCALDigit * tower ;
275 Int_t lastIndex = lastPreShowerIndex +1 ;
278 for (index = firstPreShowerIndex ; index <= lastPreShowerIndex; index++) {
279 preshower = dynamic_cast<AliEMCALDigit *>(sdigits->At(index) );
280 Bool_t towerFound = kFALSE ;
282 for (jndex = 0; jndex < firstPreShowerIndex; jndex++) {
283 tower = dynamic_cast<AliEMCALDigit *>(sdigits->At(jndex) );
284 if ( (preshower->GetId() - (geom->GetNZ() * geom->GetNPhi()) ) == tower->GetId() ) {
285 Float_t towerEnergy = static_cast<Float_t>(tower->GetAmp()) ;
286 Float_t preshoEnergy = static_cast<Float_t>(preshower->GetAmp()) ;
287 towerEnergy +=preshoEnergy ;
288 *tower = *tower + *preshower ; // and add preshower multiplied by layer ratio to tower
289 tower->SetAmp(static_cast<Int_t>(TMath::Ceil(towerEnergy))) ;
295 new((*sdigits)[lastIndex]) AliEMCALDigit(*preshower);
296 AliEMCALDigit * temp = dynamic_cast<AliEMCALDigit *>(sdigits->At(lastIndex)) ;
297 temp->SetId(temp->GetId() - (geom->GetNZ() * geom->GetNPhi()) ) ;
303 Int_t NPrimarymax = -1 ;
304 for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) {
305 sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(i)) ;
306 sdigit->SetIndexInList(i) ;
309 for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) {
310 if (((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) > NPrimarymax)
311 NPrimarymax = ((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) ;
314 if(gAlice->TreeS() == 0)
315 gAlice->MakeTree("S",fSplitFile);
317 //First list of sdigits
318 Int_t bufferSize = 32000 ;
319 TBranch * sdigitsBranch = gAlice->TreeS()->Branch("EMCAL",&sdigits,bufferSize);
320 sdigitsBranch->SetTitle(sdname);
322 //second - SDigitizer
323 Int_t splitlevel = 0 ;
324 AliEMCALSDigitizer * sd = this ;
325 TBranch * sdigitizerBranch = gAlice->TreeS()->Branch("AliEMCALSDigitizer","AliEMCALSDigitizer",
326 &sd,bufferSize,splitlevel);
327 sdigitizerBranch->SetTitle(sdname);
329 sdigitsBranch->Fill() ;
330 sdigitizerBranch->Fill() ;
331 gAlice->TreeS()->AutoSave() ;
333 if(strstr(option,"deb"))
334 PrintSDigits(option) ;
339 if ( fSplitFile->IsOpen() )
340 fSplitFile->Close() ;
342 if(strstr(option,"tim")){
343 gBenchmark->Stop("EMCALSDigitizer");
344 cout << "AliEMCALSDigitizer:" << endl ;
345 cout << " took " << gBenchmark->GetCpuTime("EMCALSDigitizer") << " seconds for SDigitizing "
346 << gBenchmark->GetCpuTime("EMCALSDigitizer")/fNevents << " seconds per event " << endl ;
352 //__________________________________________________________________
353 void AliEMCALSDigitizer::SetSDigitsBranch(const char * title ){
355 // Setting title to branch SDigits
357 TString stitle(title) ;
359 // check if branch with title already exists
360 TBranch * sdigitsBranch =
361 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("EMCAL")) ;
362 TBranch * sdigitizerBranch =
363 static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("AliEMCALSDigitizer")) ;
364 const char * sdigitsTitle = sdigitsBranch ->GetTitle() ;
365 const char * sdigitizerTitle = sdigitizerBranch ->GetTitle() ;
366 if ( stitle.CompareTo(sdigitsTitle)==0 || stitle.CompareTo(sdigitizerTitle)==0 ){
367 cerr << "ERROR: AliEMCALSdigitizer::SetSDigitsBranch -> Cannot overwrite existing branch with title " << title << endl ;
371 cout << "AliEMCALSdigitizer::SetSDigitsBranch -> Changing SDigits file from " << GetName() << " to " << title << endl ;
375 // Post to the WhiteBoard
376 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
377 gime->PostSDigits( title, GetTitle()) ;
382 //__________________________________________________________________
383 void AliEMCALSDigitizer::SetSplitFile(const TString splitFileName)
385 // Diverts the SDigits in a file separate from the hits file
387 TDirectory * cwd = gDirectory ;
389 if ( !(gAlice->GetTreeSFileName() == splitFileName) ) {
390 if (gAlice->GetTreeSFile() )
391 gAlice->GetTreeSFile()->Close() ;
394 fSplitFile = gAlice->InitTreeFile("S",splitFileName.Data());
396 if ( !fSplitFile->Get("gAlice") )
399 TTree *treeE = gAlice->TreeE();
401 cerr << "ERROR: AliEMCALSDigitizer::SetSPlitFile -> No TreeE found "<<endl;
406 if ( !fSplitFile->Get("TreeE") ) {
407 AliHeader *header = new AliHeader();
408 treeE->SetBranchAddress("Header", &header);
409 treeE->SetBranchStatus("*",1);
410 TTree *treeENew = treeE->CloneTree();
415 if ( !fSplitFile->Get("AliceGeom") ) {
416 TGeometry *AliceGeom = static_cast<TGeometry*>(cwd->Get("AliceGeom"));
418 cerr << "ERROR: AliEMCALSDigitizer::SetSPlitFile -> AliceGeom was not found in the input file "<<endl;
424 gAlice->MakeTree("S",fSplitFile);
426 cout << "INFO: AliEMCALSDigitizer::SetSPlitFile -> SDigits will be stored in " << splitFileName.Data() << endl ;
429 //__________________________________________________________________
430 void AliEMCALSDigitizer::Print(Option_t* option)const{
431 cout << "------------------- "<< GetName() << " -------------" << endl ;
432 cout << " Writing SDigitis to branch with title " << GetName() << endl ;
433 cout << " with digitization parameters A = " << fA << endl ;
434 cout << " B = " << fB << endl ;
435 cout << " Threshold for Primary assignment in Tower = " << fTowerPrimThreshold << endl ;
436 cout << " Threshold for Primary assignment in PreShower = " << fPreShowerPrimThreshold << endl ;
437 cout << "---------------------------------------------------"<<endl ;
440 //__________________________________________________________________
441 Bool_t AliEMCALSDigitizer::operator==( AliEMCALSDigitizer const &sd )const{
442 if( (fA==sd.fA)&&(fB==sd.fB)&&
443 (fTowerPrimThreshold==sd.fTowerPrimThreshold) &&
444 (fPreShowerPrimThreshold==sd.fPreShowerPrimThreshold))
449 //__________________________________________________________________
450 void AliEMCALSDigitizer::PrintSDigits(Option_t * option){
451 //Prints list of digits produced at the current pass of AliEMCALDigitizer
453 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
454 TString sdname(GetName()) ;
455 sdname.Remove(sdname.Index(GetTitle())-1) ;
456 const TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
458 cout << "AliEMCALSDigitiser: event " << gAlice->GetEvNumber() << endl ;
459 cout << " Number of entries in SDigits list " << sdigits->GetEntriesFast() << endl ;
461 if(strstr(option,"all")||strstr(option,"EMC")){
464 AliEMCALDigit * digit;
465 cout << "SDigit Id " << " Amplitude " << " Time " << " Index " << " Nprim " << " Primaries list " << endl;
467 for (index = 0 ; index < sdigits->GetEntries() ; index++) {
468 digit = (AliEMCALDigit * ) sdigits->At(index) ;
469 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " " << digit->GetTime()
470 << setw(6) << digit->GetIndexInList() << " "
471 << setw(5) << digit->GetNprimary() <<" ";
474 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
475 cout << " " << digit->GetPrimary(iprimary+1) << " ";
481 //________________________________________________________________________
482 Int_t AliEMCALSDigitizer::Layer2TowerID(Int_t ihit, Bool_t preshower){
483 // Method to Transform from Hit Id to Digit Id
484 // This function should be one to one
485 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
486 const AliEMCALGeometry * geom = gime->EMCALGeometry();
487 Int_t ieta = ((ihit-1)/geom->GetNPhi())%geom->GetNZ(); // eta Tower Index
488 Int_t iphi = (ihit-1)%(geom->GetNPhi())+1; //phi Tower Index
492 if (preshower)ipre = 1;
493 if (iphi > 0 && ieta >= 0){
494 it = iphi+ieta*geom->GetNPhi() + ipre*geom->GetNPhi()*geom->GetNZ();
497 cerr << " AliEMCALSDigitizer::Layer2TowerID() -- there is an error "<< endl << "Eta number = "
498 << ieta << "Phi number = " << iphi << endl ;
500 } // end if iphi>0 && ieta>0
502 //_______________________________________________________________________________________
503 // void AliEMCALSDigitizer::TestTowerID(void)
507 // Bool_t preshower = kFALSE;
508 // for (j = 0 ; j < 10 ; j++){ // loop over hit id
510 // for (i = 0 ; i <= 2 ; i++){ // loop over
511 // Int_t k = i*96*144+j*144+1;
512 // cout << " Hit Index = " << k << " " << j*10 << " TOWERID = " << Layer2TowerID(k, preshower) << endl ;