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 primary in one active cell
22 // A threshold for assignment of the primary to SDigit is applied
23 // SDigits are written to TreeS, branch "PHOS"
24 // AliPHOSSDigitizer with all current parameters is written
25 // to TreeS branch "AliPHOSSDigitizer".
26 // Both branches, "PHOS" and "AliPHOSSDigitizer", are written to the same
27 // file, and therefore, changing branch file name one can produce several
28 // versions of SDigitization from the same hits.
31 //*-- Author : Dmitri Peressounko (SUBATECH & KI)
32 //////////////////////////////////////////////////////////////////////////////
34 // --- ROOT system ---
41 #include "TBenchmark.h"
42 // --- Standard library ---
45 // --- AliRoot header files ---
47 #include "AliPHOSDigit.h"
48 #include "AliPHOSHit.h"
49 #include "AliPHOSSDigitizer.h"
52 ClassImp(AliPHOSSDigitizer)
55 //____________________________________________________________________________
56 AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("AliPHOSSDigitizer","")
61 fPrimThreshold = 0.01 ;
65 fIsInitialized = kFALSE ;
69 //____________________________________________________________________________
70 AliPHOSSDigitizer::AliPHOSSDigitizer(const char* headerFile, const char *sDigitsTitle):TTask("AliPHOSSDigitizer","")
75 fPrimThreshold = 0.01 ;
77 fSDigitsTitle = sDigitsTitle ;
78 fHeadersFile = headerFile ;
79 fSDigits = new TClonesArray("AliPHOSDigit",1000);
80 fHits = new TClonesArray("AliPHOSHit",1000);
82 TFile * file = (TFile*) gROOT->GetFile(fHeadersFile.Data() ) ;
84 //File was not opened yet
86 file = new TFile(fHeadersFile.Data(),"update") ;
87 gAlice = (AliRun *) file->Get("gAlice") ;
90 //add Task to //root/Tasks folder
91 TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ;
92 roottasks->Add(this) ;
94 fIsInitialized = kTRUE ;
97 //____________________________________________________________________________
98 AliPHOSSDigitizer::~AliPHOSSDigitizer()
106 //____________________________________________________________________________
107 void AliPHOSSDigitizer::Init(){
108 //Initialization can not be done in the default constructor
112 if(fHeadersFile.IsNull())
113 fHeadersFile="galice.root" ;
115 TFile * file = (TFile*) gROOT->GetFile(fHeadersFile.Data() ) ;
117 //if file was not opened yet, read gAlice
119 file = new TFile(fHeadersFile.Data(),"update") ;
120 gAlice = (AliRun *) file->Get("gAlice") ;
123 fHits = new TClonesArray("AliPHOSHit",1000);
124 fSDigits = new TClonesArray("AliPHOSDigit",1000);
126 // add Task to //root/Tasks folder
127 TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ;
128 roottasks->Add(this) ;
130 fIsInitialized = kTRUE ;
133 //____________________________________________________________________________
134 void AliPHOSSDigitizer::Exec(Option_t *option) {
135 //Collects all hits in the same active volume into digit
140 if(strstr(option,"tim"))
141 gBenchmark->Start("PHOSSDigitizer");
143 fNevents = (Int_t) gAlice->TreeE()->GetEntries() ;
146 for(ievent = 0; ievent < fNevents; ievent++){
147 gAlice->GetEvent(ievent) ;
148 gAlice->SetEvent(ievent) ;
150 if(gAlice->TreeH()==0){
151 cout << "AliPHOSSDigitizer: There is no Hit Tree" << endl;
155 //set address of the hits
156 TBranch * branch = gAlice->TreeH()->GetBranch("PHOS");
158 branch->SetAddress(&fHits);
160 cout << "ERROR in AliPHOSSDigitizer: "<< endl ;
161 cout << " no branch PHOS in TreeH"<< endl ;
162 cout << " do nothing " << endl ;
170 //Now made SDigits from hits, for PHOS it is the same, so just copy
172 for (itrack=0; itrack<gAlice->GetNtrack(); itrack++){
174 //=========== Get the Hits Tree for the Primary track itrack
176 gAlice->TreeH()->GetEvent(itrack);
179 for ( i = 0 ; i < fHits->GetEntries() ; i++ ) {
180 AliPHOSHit * hit = (AliPHOSHit*)fHits->At(i) ;
181 AliPHOSDigit * newdigit ;
183 // Assign primary number only if contribution is significant
184 if( hit->GetEnergy() > fPrimThreshold)
185 newdigit = new AliPHOSDigit( hit->GetPrimary(), hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
187 newdigit = new AliPHOSDigit( -1 , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
189 new((*fSDigits)[nSdigits]) AliPHOSDigit(* newdigit) ;
195 } // loop over tracks
199 nSdigits = fSDigits->GetEntriesFast() ;
200 fSDigits->Expand(nSdigits) ;
203 for (i = 0 ; i < nSdigits ; i++) {
204 AliPHOSDigit * digit = (AliPHOSDigit *) fSDigits->At(i) ;
205 digit->SetIndexInList(i) ;
208 if(gAlice->TreeS() == 0)
209 gAlice->MakeTree("S") ;
211 //check, if this branch already exits?
212 TBranch * sdigitsBranch = 0;
213 TBranch * sdigitizerBranch = 0;
215 TObjArray * branches = gAlice->TreeS()->GetListOfBranches() ;
217 Bool_t phosNotFound = kTRUE ;
218 Bool_t sdigitizerNotFound = kTRUE ;
220 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
223 sdigitsBranch=(TBranch *) branches->At(ibranch) ;
224 if( (strcmp("PHOS",sdigitsBranch->GetName())==0 ) &&
225 (fSDigitsTitle.CompareTo(sdigitsBranch->GetTitle()) == 0) )
226 phosNotFound = kFALSE ;
228 if(sdigitizerNotFound){
229 sdigitizerBranch = (TBranch *) branches->At(ibranch) ;
230 if( (strcmp(sdigitizerBranch->GetName(),"AliPHOSSDigitizer") == 0)&&
231 (fSDigitsTitle.CompareTo(sdigitizerBranch->GetTitle()) == 0) )
232 sdigitizerNotFound = kFALSE ;
236 if(!(sdigitizerNotFound && phosNotFound)){
237 cout << "AliPHOSSdigitizer error:" << endl ;
238 cout << "Can not overwrite existing branches: do not write" << endl ;
242 //Make (if necessary) branches
244 if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
245 file = new char[strlen(gAlice->GetBaseFile())+20] ;
246 sprintf(file,"%s/PHOS.SDigits.root",gAlice->GetBaseFile()) ;
249 TDirectory *cwd = gDirectory;
251 //First list of sdigits
252 Int_t bufferSize = 32000 ;
253 sdigitsBranch = gAlice->TreeS()->Branch("PHOS",&fSDigits,bufferSize);
254 sdigitsBranch->SetTitle(fSDigitsTitle.Data());
256 sdigitsBranch->SetFile(file);
257 TIter next( sdigitsBranch->GetListOfBranches());
258 while ((sdigitsBranch=(TBranch*)next())) {
259 sdigitsBranch->SetFile(file);
264 //second - SDigitizer
265 Int_t splitlevel = 0 ;
266 AliPHOSSDigitizer * sd = this ;
267 sdigitizerBranch = gAlice->TreeS()->Branch("AliPHOSSDigitizer","AliPHOSSDigitizer",
268 &sd,bufferSize,splitlevel);
269 sdigitizerBranch->SetTitle(fSDigitsTitle.Data());
271 sdigitizerBranch->SetFile(file);
272 TIter next( sdigitizerBranch->GetListOfBranches());
273 while ((sdigitizerBranch=(TBranch*)next())) {
274 sdigitizerBranch->SetFile(file);
280 gAlice->TreeS()->Fill() ;
281 gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
283 if(strstr(option,"deb"))
284 PrintSDigits(option) ;
288 if(strstr(option,"tim")){
289 gBenchmark->Stop("PHOSSDigitizer");
290 cout << "AliPHOSSDigitizer:" << endl ;
291 cout << " took " << gBenchmark->GetCpuTime("PHOSSDigitizer") << " seconds for SDigitizing "
292 << gBenchmark->GetCpuTime("PHOSSDigitizer")/fNevents << " seconds per event " << endl ;
298 //__________________________________________________________________
299 void AliPHOSSDigitizer::SetSDigitsBranch(const char * title ){
300 //Seting title to branch SDigits
301 if(!fSDigitsTitle.IsNull())
302 cout << "AliPHOSSdigitizer: changing SDigits file from " <<fSDigitsTitle.Data() << " to " << title << endl ;
303 fSDigitsTitle=title ;
305 //__________________________________________________________________
306 void AliPHOSSDigitizer::Print(Option_t* option)const{
307 cout << "------------------- "<< GetName() << " -------------" << endl ;
308 cout << " Writing SDigitis to branch with title " << fSDigitsTitle.Data() << endl ;
309 cout << " with digitization parameters A = " << fA << endl ;
310 cout << " B = " << fB << endl ;
311 cout << " Threshold for Primary assignment= " << fPrimThreshold << endl ;
312 cout << "---------------------------------------------------"<<endl ;
315 //__________________________________________________________________
316 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const{
317 if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
322 //__________________________________________________________________
323 void AliPHOSSDigitizer::PrintSDigits(Option_t * option){
324 //Prints list of digits produced at the current pass of AliPHOSDigitizer
326 cout << "AliPHOSSDigitizer: " << endl ;
327 cout << " Number of entries in SDigits list " << fSDigits->GetEntriesFast() << endl ;
330 if(strstr(option,"all")){// print all digits
333 AliPHOSDigit * digit;
334 cout << "SDigit Id " << " Amplitude " << " Index " << " Nprim " << " Primaries list " << endl;
336 for (index = 0 ; index < fSDigits->GetEntries() ; index++) {
337 digit = (AliPHOSDigit * ) fSDigits->At(index) ;
338 cout << setw(8) << digit->GetId() << " " << setw(3) << digit->GetAmp() << " "
339 << setw(6) << digit->GetIndexInList() << " "
340 << setw(5) << digit->GetNprimary() <<" ";
343 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
344 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";