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
47 // August 2002 Yves Schutz: clone PHOS as closely as possible and intoduction
48 // of new IO (à la PHOS)
49 //////////////////////////////////////////////////////////////////////////////
51 // --- ROOT system ---
58 #include "TBenchmark.h"
59 #include "TGeometry.h"
61 // --- Standard library ---
62 #include <Riostream.h>
64 // --- AliRoot header files ---
66 #include "AliHeader.h"
67 #include "AliEMCALDigit.h"
68 #include "AliEMCALGeometry.h"
69 #include "AliEMCALGetter.h"
70 #include "AliEMCALHit.h"
71 #include "AliEMCALSDigitizer.h"
73 ClassImp(AliEMCALSDigitizer)
75 //____________________________________________________________________________
76 AliEMCALSDigitizer::AliEMCALSDigitizer():TTask("AliEMCALSDigitizer","")
80 fDefaultInit = kTRUE ;
83 //____________________________________________________________________________
84 AliEMCALSDigitizer::AliEMCALSDigitizer(const char* headerFile, const char *sDigitsTitle, const Bool_t toSplit):
85 TTask(sDigitsTitle, headerFile)
91 fDefaultInit = kFALSE ;
94 //____________________________________________________________________________
95 AliEMCALSDigitizer::~AliEMCALSDigitizer()
101 //____________________________________________________________________________
102 void AliEMCALSDigitizer::Init(){
103 // Initialization: open root-file, allocate arrays for hits and sdigits,
104 // attach task SDigitizer to the list of EMCAL tasks
106 // Initialization can not be done in the default constructor
107 //============================================================= YS
108 // The initialisation is now done by the getter
110 if( strcmp(GetTitle(), "") == 0 )
111 SetTitle("galice.root") ;
113 AliEMCALGetter * gime = AliEMCALGetter::GetInstance(GetTitle(), GetName(), fToSplit) ;
115 cerr << "ERROR: AliEMCALSDigitizer::Init -> Could not obtain the Getter object !"
120 gime->PostSDigits( GetName(), GetTitle() ) ;
124 // construct the name of the file as /path/EMCAL.SDigits.root
125 // First - extract full path if necessary
126 TString sDigitsFileName(GetTitle()) ;
127 Ssiz_t islash = sDigitsFileName.Last('/') ;
128 if(islash<sDigitsFileName.Length())
129 sDigitsFileName.Remove(islash+1,sDigitsFileName.Length()) ;
133 // Next - append the file name
134 sDigitsFileName+="EMCAL.SDigits." ;
135 if((strcmp(GetName(),"Default")!=0)&&(strcmp(GetName(),"")!=0)){
136 sDigitsFileName+=GetName() ;
137 sDigitsFileName+="." ;
139 sDigitsFileName+="root" ;
141 // Finally - check if the file already opened or open the file
142 fSplitFile = static_cast<TFile*>(gROOT->GetFile(sDigitsFileName.Data()));
144 fSplitFile = TFile::Open(sDigitsFileName.Data(),"update") ;
147 TString sdname(GetName() );
149 sdname.Append(GetTitle() ) ;
151 gime->PostSDigitizer(this) ;
154 //____________________________________________________________________________
155 void AliEMCALSDigitizer::InitParameters()
159 fTowerPrimThreshold = 0.01 ;
160 fPreShowerPrimThreshold = 0.0001 ;
161 fPhotonElectronFactor = 5000. ; // photoelectrons per GeV
166 //____________________________________________________________________________
167 void AliEMCALSDigitizer::Exec(Option_t *option)
169 // Collects all hits in the same active volume into digit
171 if( strcmp(GetName(), "") == 0 )
174 if (strstr(option, "print") ) {
179 if(strstr(option,"tim"))
180 gBenchmark->Start("EMCALSDigitizer");
182 //Check, if this branch already exits
183 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
184 if(gime->BranchExists("SDigits") )
187 TString sdname(GetName()) ;
188 sdname.Remove(sdname.Index(GetTitle())-1) ;
190 Int_t nevents = gime->MaxEvent() ;
193 for(ievent = 0; ievent < nevents; ievent++){
194 gime->Event(ievent,"H") ;
196 const TClonesArray * hits = gime->Hits() ;
198 TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
202 //Collects all hits in the same active volume into digit
204 //Now make SDigits from hits, for EMCAL it is the same, so just copy
205 Int_t nPrim = static_cast<Int_t>((gAlice->TreeH())->GetEntries()) ;
206 // Attention nPrim is the number of primaries tracked by Geant
207 // and this number could be different to the number of Primaries in TreeK;
209 for ( iprim = 0 ; iprim < nPrim ; iprim++ ) {
210 //=========== Get the EMCAL branch from Hits Tree for the Primary iprim
213 for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
214 AliEMCALHit * hit = dynamic_cast<AliEMCALHit*>(hits->At(i)) ;
215 AliEMCALDigit * curSDigit = 0 ;
216 AliEMCALDigit * sdigit = 0 ;
217 Bool_t newsdigit = kTRUE;
218 // Assign primary number only if deposited energy is significant
220 if( (!hit->IsInPreShower() && hit->GetEnergy() > fTowerPrimThreshold) ||
221 (hit->IsInPreShower() && hit->GetEnergy() > fPreShowerPrimThreshold))
222 curSDigit = new AliEMCALDigit( hit->GetPrimary(),
223 hit->GetIparent(),Layer2TowerID(hit->GetId(),hit->IsInPreShower()),
224 Digitize(hit->GetEnergy()), hit->GetTime() ) ;
226 curSDigit = new AliEMCALDigit( -1 ,
228 Layer2TowerID(hit->GetId(),hit->IsInPreShower()),
229 Digitize(hit->GetEnergy()), hit->GetTime() ) ;
231 for(check= 0; check < nSdigits ; check++) {
232 sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(check)) ;
233 if( sdigit->GetId() == curSDigit->GetId()) { // Are we in the same tower or the same preshower ?
234 *sdigit = *sdigit + *curSDigit;
239 new((*sdigits)[nSdigits]) AliEMCALDigit(*curSDigit);
243 } // loop over all hits (hit = deposited energy/layer/entering particle)
248 nSdigits = sdigits->GetEntriesFast() ;
249 fSDigitsInRun += nSdigits ;
250 sdigits->Expand(nSdigits) ;
252 const AliEMCALGeometry * geom = gime->EMCALGeometry() ;
254 if (nSdigits != 0 ) {
255 Int_t lastPreShowerIndex = nSdigits - 1 ;
258 if (!(dynamic_cast<AliEMCALDigit *>(sdigits->At(lastPreShowerIndex))->IsInPreShower()))
260 lastPreShowerIndex = -2;
262 Int_t firstPreShowerIndex = 100000 ;
264 AliEMCALDigit * sdigit = 0 ;
266 for ( index = 0; index < nSdigits ; index++) {
267 sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(index) ) ;
268 if (sdigit->IsInPreShower() ){
269 firstPreShowerIndex = index ;
274 AliEMCALDigit * preshower ;
275 AliEMCALDigit * tower ;
276 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))) ;
294 new((*sdigits)[lastIndex]) AliEMCALDigit(*preshower);
295 AliEMCALDigit * temp = dynamic_cast<AliEMCALDigit *>(sdigits->At(lastIndex)) ;
296 temp->SetId(temp->GetId() - (geom->GetNZ() * geom->GetNPhi()) ) ;
301 Int_t NPrimarymax = -1 ;
303 for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) {
304 sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(i)) ;
305 sdigit->SetIndexInList(i) ;
308 for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) {
309 if (((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) > NPrimarymax)
310 NPrimarymax = ((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) ;
316 if(gAlice->TreeS() == 0 || (fSplitFile)) //<--- To be checked: we should not create TreeS if it is already here
317 gAlice->MakeTree("S",fSplitFile);
322 //First list of sdigits
323 Int_t bufferSize = 32000 ;
324 TBranch * sdigitsBranch = gAlice->TreeS()->Branch("EMCAL",&sdigits,bufferSize);
325 sdigitsBranch->SetTitle(sdname);
328 Int_t splitlevel = 0 ;
329 AliEMCALSDigitizer * sd = this ;
330 TBranch * sdigitizerBranch = gAlice->TreeS()->Branch("AliEMCALSDigitizer","AliEMCALSDigitizer",
331 &sd,bufferSize,splitlevel);
332 sdigitizerBranch->SetTitle(sdname);
334 sdigitsBranch->Fill() ;
335 sdigitizerBranch->Fill() ;
336 gAlice->TreeS()->AutoSave() ;
338 if(strstr(option,"deb"))
339 PrintSDigits(option) ;
343 if(strstr(option,"tim")){
344 gBenchmark->Stop("EMCALSDigitizer");
345 cout << "AliEMCALSDigitizer:" << endl ;
346 cout << " took " << gBenchmark->GetCpuTime("EMCALSDigitizer") << " seconds for SDigitizing "
347 << gBenchmark->GetCpuTime("EMCALSDigitizer") << " 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()) ;
381 //__________________________________________________________________
382 void AliEMCALSDigitizer::Print(Option_t* option)const
384 // Prints parameters of SDigitizer
386 cout << "------------------- "<< GetName() << " -------------" << endl ;
387 cout << " Writing SDigitis to branch with title " << GetName() << endl ;
388 cout << " with digitization parameters A = " << fA << endl ;
389 cout << " B = " << fB << endl ;
390 cout << " Threshold for Primary assignment in Tower = " << fTowerPrimThreshold << endl ;
391 cout << " Threshold for Primary assignment in PreShower = " << fPreShowerPrimThreshold << endl ;
392 cout << "---------------------------------------------------"<<endl ;
396 //__________________________________________________________________
397 Bool_t AliEMCALSDigitizer::operator==( AliEMCALSDigitizer const &sd )const
400 // SDititizers are equal if their pedestal, slope and threshold are equal
402 if( (fA==sd.fA)&&(fB==sd.fB)&&
403 (fTowerPrimThreshold==sd.fTowerPrimThreshold) &&
404 (fPreShowerPrimThreshold==sd.fPreShowerPrimThreshold))
410 //__________________________________________________________________
411 void AliEMCALSDigitizer::PrintSDigits(Option_t * option){
412 //Prints list of digits produced at the current pass of AliEMCALDigitizer
414 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
415 TString sdname(GetName()) ;
416 sdname.Remove(sdname.Index(GetTitle())-1) ;
417 const TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
419 cout << "AliEMCALSDigitiser: event " << gAlice->GetEvNumber() << endl ;
420 cout << " Number of entries in SDigits list " << sdigits->GetEntriesFast() << endl ;
422 if(strstr(option,"all")||strstr(option,"EMC")){
425 AliEMCALDigit * digit;
426 cout << "SDigit Id " << " Amplitude " << " Time " << " Index " << " Nprim " << " Primaries list " << endl;
428 for (index = 0 ; index < sdigits->GetEntries() ; index++) {
429 digit = dynamic_cast<AliEMCALDigit *>( sdigits->At(index) ) ;
430 cout << setw(6) << digit->GetId() << " " << setw(10) << digit->GetAmp() << " " << digit->GetTime()
431 << setw(6) << digit->GetIndexInList() << " "
432 << setw(5) << digit->GetNprimary() <<" ";
435 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
436 cout << " " << digit->GetPrimary(iprimary+1) << " ";
443 //________________________________________________________________________
444 const Int_t AliEMCALSDigitizer::Layer2TowerID(Int_t ihit, Bool_t preshower)
446 // Method to Transform from Hit Id to Digit Id
447 // This function should be one to one
448 AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ;
449 const AliEMCALGeometry * geom = gime->EMCALGeometry();
450 Int_t ieta = ((ihit-1)/geom->GetNPhi())%geom->GetNZ(); // eta Tower Index
451 Int_t iphi = (ihit-1)%(geom->GetNPhi())+1; //phi Tower Index
455 if (preshower)ipre = 1;
456 if (iphi > 0 && ieta >= 0){
457 it = iphi+ieta*geom->GetNPhi() + ipre*geom->GetNPhi()*geom->GetNZ();
460 cerr << " AliEMCALSDigitizer::Layer2TowerID() -- there is an error "<< endl << "Eta number = "
461 << ieta << "Phi number = " << iphi << endl ;
463 } // end if iphi>0 && ieta>0
465 //_______________________________________________________________________________________
466 // void AliEMCALSDigitizer::TestTowerID(void)
470 // Bool_t preshower = kFALSE;
471 // for (j = 0 ; j < 10 ; j++){ // loop over hit id
473 // for (i = 0 ; i <= 2 ; i++){ // loop over
474 // Int_t k = i*96*144+j*144+1;
475 // cout << " Hit Index = " << k << " " << j*10 << " TOWERID = " << Layer2TowerID(k, preshower) << endl ;
480 //____________________________________________________________________________
481 void AliEMCALSDigitizer::UseHitsFrom(const char * filename)