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 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] AliPHOSSDigitizer * s = new AliPHOSSDigitizer("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 : Dmitri Peressounko (SUBATECH & KI)
45 //////////////////////////////////////////////////////////////////////////////
48 // --- ROOT system ---
55 #include "TBenchmark.h"
56 // --- Standard library ---
59 // --- AliRoot header files ---
61 #include "AliPHOSDigit.h"
62 #include "AliPHOSHit.h"
63 #include "AliPHOSSDigitizer.h"
66 ClassImp(AliPHOSSDigitizer)
69 //____________________________________________________________________________
70 AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("AliPHOSSDigitizer","")
75 fPrimThreshold = 0.01 ;
79 fIsInitialized = kFALSE ;
83 //____________________________________________________________________________
84 AliPHOSSDigitizer::AliPHOSSDigitizer(const char* headerFile, const char *sDigitsTitle):TTask("AliPHOSSDigitizer","")
89 fPrimThreshold = 0.01 ;
91 fSDigitsTitle = sDigitsTitle ;
92 fHeadersFile = headerFile ;
93 fSDigits = new TClonesArray("AliPHOSDigit",1000);
94 fHits = new TClonesArray("AliPHOSHit",1000);
96 TFile * file = (TFile*) gROOT->GetFile(fHeadersFile.Data() ) ;
98 //File was not opened yet
100 if(fHeadersFile.Contains("rfio"))
101 file = TFile::Open(fHeadersFile,"update") ;
103 file = new TFile(fHeadersFile.Data(),"update") ;
104 gAlice = (AliRun *) file->Get("gAlice") ;
107 //add Task to //root/Tasks folder
108 TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ;
109 roottasks->Add(this) ;
111 fIsInitialized = kTRUE ;
114 //____________________________________________________________________________
115 AliPHOSSDigitizer::~AliPHOSSDigitizer()
123 //____________________________________________________________________________
124 void AliPHOSSDigitizer::Init(){
125 //Initialization can not be done in the default constructor
129 if(fHeadersFile.IsNull())
130 fHeadersFile="galice.root" ;
132 TFile * file = (TFile*) gROOT->GetFile(fHeadersFile.Data() ) ;
134 //if file was not opened yet, read gAlice
136 file = new TFile(fHeadersFile.Data(),"update") ;
137 gAlice = (AliRun *) file->Get("gAlice") ;
140 fHits = new TClonesArray("AliPHOSHit",1000);
141 fSDigits = new TClonesArray("AliPHOSDigit",1000);
143 // add Task to //root/Tasks folder
144 TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ;
145 roottasks->Add(this) ;
147 fIsInitialized = kTRUE ;
150 //____________________________________________________________________________
151 void AliPHOSSDigitizer::Exec(Option_t *option) {
152 //Collects all hits in the same active volume into digit
157 if(strstr(option,"tim"))
158 gBenchmark->Start("PHOSSDigitizer");
160 fNevents = (Int_t) gAlice->TreeE()->GetEntries() ;
163 for(ievent = 0; ievent < fNevents; ievent++){
164 gAlice->GetEvent(ievent) ;
165 gAlice->SetEvent(ievent) ;
167 if(gAlice->TreeH()==0){
168 cout << "AliPHOSSDigitizer: There is no Hit Tree" << endl;
172 //set address of the hits
173 TBranch * branch = gAlice->TreeH()->GetBranch("PHOS");
175 branch->SetAddress(&fHits);
177 cout << "ERROR in AliPHOSSDigitizer: "<< endl ;
178 cout << " no branch PHOS in TreeH"<< endl ;
179 cout << " do nothing " << endl ;
187 //Now made SDigits from hits, for PHOS it is the same, so just copy
189 for (itrack=0; itrack < gAlice->GetNtrack(); itrack++){
191 //=========== Get the PHOS branch from Hits Tree for the Primary track itrack
192 branch->GetEntry(itrack,0);
195 for ( i = 0 ; i < fHits->GetEntries() ; i++ ) {
196 AliPHOSHit * hit = (AliPHOSHit*)fHits->At(i) ;
198 // Assign primary number only if contribution is significant
199 if( hit->GetEnergy() > fPrimThreshold)
200 new((*fSDigits)[nSdigits]) AliPHOSDigit( hit->GetPrimary(), hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
202 new((*fSDigits)[nSdigits]) AliPHOSDigit( -1 , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
208 } // loop over tracks
212 nSdigits = fSDigits->GetEntriesFast() ;
213 fSDigits->Expand(nSdigits) ;
216 for (i = 0 ; i < nSdigits ; i++) {
217 AliPHOSDigit * digit = (AliPHOSDigit *) fSDigits->At(i) ;
218 digit->SetIndexInList(i) ;
221 if(gAlice->TreeS() == 0)
222 gAlice->MakeTree("S") ;
224 //check, if this branch already exits?
225 TBranch * sdigitsBranch = 0;
226 TBranch * sdigitizerBranch = 0;
228 TObjArray * branches = gAlice->TreeS()->GetListOfBranches() ;
230 Bool_t phosNotFound = kTRUE ;
231 Bool_t sdigitizerNotFound = kTRUE ;
233 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
236 sdigitsBranch=(TBranch *) branches->At(ibranch) ;
237 if( (strcmp("PHOS",sdigitsBranch->GetName())==0 ) &&
238 (fSDigitsTitle.CompareTo(sdigitsBranch->GetTitle()) == 0) )
239 phosNotFound = kFALSE ;
241 if(sdigitizerNotFound){
242 sdigitizerBranch = (TBranch *) branches->At(ibranch) ;
243 if( (strcmp(sdigitizerBranch->GetName(),"AliPHOSSDigitizer") == 0)&&
244 (fSDigitsTitle.CompareTo(sdigitizerBranch->GetTitle()) == 0) )
245 sdigitizerNotFound = kFALSE ;
249 if(!(sdigitizerNotFound && phosNotFound)){
250 cout << "AliPHOSSdigitizer error:" << endl ;
251 cout << "Can not overwrite existing branches: do not write" << endl ;
255 //Make (if necessary) branches
257 if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
258 file = new char[strlen(gAlice->GetBaseFile())+20] ;
259 sprintf(file,"%s/PHOS.SDigits.root",gAlice->GetBaseFile()) ;
262 TDirectory *cwd = gDirectory;
264 //First list of sdigits
265 Int_t bufferSize = 32000 ;
266 sdigitsBranch = gAlice->TreeS()->Branch("PHOS",&fSDigits,bufferSize);
267 sdigitsBranch->SetTitle(fSDigitsTitle.Data());
269 sdigitsBranch->SetFile(file);
270 TIter next( sdigitsBranch->GetListOfBranches());
271 while ((sdigitsBranch=(TBranch*)next())) {
272 sdigitsBranch->SetFile(file);
277 //second - SDigitizer
278 Int_t splitlevel = 0 ;
279 AliPHOSSDigitizer * sd = this ;
280 sdigitizerBranch = gAlice->TreeS()->Branch("AliPHOSSDigitizer","AliPHOSSDigitizer",
281 &sd,bufferSize,splitlevel);
282 sdigitizerBranch->SetTitle(fSDigitsTitle.Data());
284 sdigitizerBranch->SetFile(file);
285 TIter next( sdigitizerBranch->GetListOfBranches());
286 while ((sdigitizerBranch=(TBranch*)next())) {
287 sdigitizerBranch->SetFile(file);
293 gAlice->TreeS()->Fill() ;
294 gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
296 if(strstr(option,"deb"))
297 PrintSDigits(option) ;
301 if(strstr(option,"tim")){
302 gBenchmark->Stop("PHOSSDigitizer");
303 cout << "AliPHOSSDigitizer:" << endl ;
304 cout << " took " << gBenchmark->GetCpuTime("PHOSSDigitizer") << " seconds for SDigitizing "
305 << gBenchmark->GetCpuTime("PHOSSDigitizer")/fNevents << " seconds per event " << endl ;
311 //__________________________________________________________________
312 void AliPHOSSDigitizer::SetSDigitsBranch(const char * title ){
313 //Seting title to branch SDigits
314 if(!fSDigitsTitle.IsNull())
315 cout << "AliPHOSSdigitizer: changing SDigits file from " <<fSDigitsTitle.Data() << " to " << title << endl ;
316 fSDigitsTitle=title ;
318 //__________________________________________________________________
319 void AliPHOSSDigitizer::Print(Option_t* option)const{
320 cout << "------------------- "<< GetName() << " -------------" << endl ;
321 cout << " Writing SDigitis to branch with title " << fSDigitsTitle.Data() << endl ;
322 cout << " with digitization parameters A = " << fA << endl ;
323 cout << " B = " << fB << endl ;
324 cout << " Threshold for Primary assignment= " << fPrimThreshold << endl ;
325 cout << "---------------------------------------------------"<<endl ;
328 //__________________________________________________________________
329 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const{
330 if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
335 //__________________________________________________________________
336 void AliPHOSSDigitizer::PrintSDigits(Option_t * option){
337 //Prints list of digits produced at the current pass of AliPHOSDigitizer
339 cout << "AliPHOSSDigitizer: " << endl ;
340 cout << " Number of entries in SDigits list " << fSDigits->GetEntriesFast() << endl ;
343 if(strstr(option,"all")){// print all digits
346 AliPHOSDigit * digit;
347 cout << "SDigit Id " << " Amplitude " << " Index " << " Nprim " << " Primaries list " << endl;
349 for (index = 0 ; index < fSDigits->GetEntries() ; index++) {
350 digit = (AliPHOSDigit * ) fSDigits->At(index) ;
351 cout << setw(8) << digit->GetId() << " " << setw(3) << digit->GetAmp() << " "
352 << setw(6) << digit->GetIndexInList() << " "
353 << setw(5) << digit->GetNprimary() <<" ";
356 for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
357 cout << setw(5) << digit->GetPrimary(iprimary+1) << " ";