]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSSDigitizer.cxx
Read not all TreeH, but only PHOS branch now -> faster 300 times
[u/mrichter/AliRoot.git] / PHOS / AliPHOSSDigitizer.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id$ */
17
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.
29 // User case:
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
43 //
44 //*-- Author :  Dmitri Peressounko (SUBATECH & KI) 
45 //////////////////////////////////////////////////////////////////////////////
46
47
48 // --- ROOT system ---
49 #include "TFile.h"
50 #include "TTask.h"
51 #include "TTree.h"
52 #include "TSystem.h"
53 #include "TROOT.h"
54 #include "TFolder.h"
55 #include "TBenchmark.h"
56 // --- Standard library ---
57 #include <iomanip.h>
58
59 // --- AliRoot header files ---
60 #include "AliRun.h"
61 #include "AliPHOSDigit.h"
62 #include "AliPHOSHit.h"
63 #include "AliPHOSSDigitizer.h"
64
65
66 ClassImp(AliPHOSSDigitizer)
67
68            
69 //____________________________________________________________________________ 
70   AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("AliPHOSSDigitizer","") 
71 {
72   // ctor
73   fA = 0;
74   fB = 10000000. ;
75   fPrimThreshold = 0.01 ;
76   fNevents = 0 ;     
77   fSDigits = 0 ;
78   fHits = 0 ;
79   fIsInitialized = kFALSE ;
80
81 }
82
83 //____________________________________________________________________________ 
84 AliPHOSSDigitizer::AliPHOSSDigitizer(const char* headerFile, const char *sDigitsTitle):TTask("AliPHOSSDigitizer","")
85 {
86   // ctor
87   fA = 0;
88   fB = 10000000.;
89   fPrimThreshold = 0.01 ;
90   fNevents = 0 ;      
91   fSDigitsTitle = sDigitsTitle ;
92   fHeadersFile = headerFile ;
93   fSDigits = new TClonesArray("AliPHOSDigit",1000);
94   fHits    = new TClonesArray("AliPHOSHit",1000);
95
96   TFile * file = (TFile*) gROOT->GetFile(fHeadersFile.Data() ) ;
97   
98   //File was not opened yet
99   if(file == 0){
100     if(fHeadersFile.Contains("rfio"))
101       file =    TFile::Open(fHeadersFile,"update") ;
102     else
103       file = new TFile(fHeadersFile.Data(),"update") ;
104     gAlice = (AliRun *) file->Get("gAlice") ;
105   }
106   
107   //add Task to //root/Tasks folder
108   TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
109   roottasks->Add(this) ; 
110   
111   fIsInitialized = kTRUE ;
112 }
113
114 //____________________________________________________________________________ 
115 AliPHOSSDigitizer::~AliPHOSSDigitizer()
116 {
117   // dtor
118   if(fSDigits)
119     delete fSDigits ;
120   if(fHits)
121     delete fHits ;
122 }
123 //____________________________________________________________________________ 
124 void AliPHOSSDigitizer::Init(){
125   //Initialization can not be done in the default constructor
126
127   if(!fIsInitialized){
128
129     if(fHeadersFile.IsNull())
130       fHeadersFile="galice.root" ;
131
132     TFile * file = (TFile*) gROOT->GetFile(fHeadersFile.Data() ) ;
133     
134     //if file was not opened yet, read gAlice
135     if(file == 0){
136       file = new TFile(fHeadersFile.Data(),"update") ;
137       gAlice = (AliRun *) file->Get("gAlice") ;
138     }
139     
140     fHits    = new TClonesArray("AliPHOSHit",1000);
141     fSDigits = new TClonesArray("AliPHOSDigit",1000);
142     
143     // add Task to //root/Tasks folder
144     TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
145     roottasks->Add(this) ; 
146     
147     fIsInitialized = kTRUE ;
148   }
149 }
150 //____________________________________________________________________________
151 void AliPHOSSDigitizer::Exec(Option_t *option) { 
152   //Collects all hits in the same active volume into digit
153   
154   if(!fIsInitialized)
155     Init() ;
156
157   if(strstr(option,"tim"))
158     gBenchmark->Start("PHOSSDigitizer");
159   
160   fNevents = (Int_t) gAlice->TreeE()->GetEntries() ; 
161   
162   Int_t ievent ;
163   for(ievent = 0; ievent < fNevents; ievent++){
164     gAlice->GetEvent(ievent) ;
165     gAlice->SetEvent(ievent) ;
166     
167     if(gAlice->TreeH()==0){
168       cout << "AliPHOSSDigitizer: There is no Hit Tree" << endl;
169       return ;
170     }
171     
172     //set address of the hits 
173     TBranch * branch = gAlice->TreeH()->GetBranch("PHOS");
174     if (branch) 
175       branch->SetAddress(&fHits);
176     else{
177       cout << "ERROR in AliPHOSSDigitizer: "<< endl ;
178       cout << "      no branch PHOS in TreeH"<< endl ;
179       cout << "      do nothing " << endl ;
180       return ;
181     }
182     
183     fSDigits->Clear();
184     Int_t nSdigits = 0 ;
185     
186     
187     //Now made SDigits from hits, for PHOS it is the same, so just copy    
188     Int_t itrack ;
189     for (itrack=0; itrack < gAlice->GetNtrack(); itrack++){
190       
191       //=========== Get the PHOS branch from Hits Tree for the Primary track itrack
192       branch->GetEntry(itrack,0);
193       
194       Int_t i;
195       for ( i = 0 ; i < fHits->GetEntries() ; i++ ) {
196         AliPHOSHit * hit = (AliPHOSHit*)fHits->At(i) ;
197
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() ) ) ;
201         else
202           new((*fSDigits)[nSdigits]) AliPHOSDigit( -1               , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
203         
204         nSdigits++ ;  
205         
206       } 
207       
208     } // loop over tracks
209     
210     fSDigits->Sort() ;
211     
212     nSdigits = fSDigits->GetEntriesFast() ;
213     fSDigits->Expand(nSdigits) ;
214     
215     Int_t i ;
216     for (i = 0 ; i < nSdigits ; i++) { 
217       AliPHOSDigit * digit = (AliPHOSDigit *) fSDigits->At(i) ; 
218       digit->SetIndexInList(i) ;     
219     }
220
221     if(gAlice->TreeS() == 0)
222       gAlice->MakeTree("S") ;
223     
224     //check, if this branch already exits?
225     TBranch * sdigitsBranch = 0;
226     TBranch * sdigitizerBranch = 0;
227     
228     TObjArray * branches = gAlice->TreeS()->GetListOfBranches() ;
229     Int_t ibranch;
230     Bool_t phosNotFound = kTRUE ;
231     Bool_t sdigitizerNotFound = kTRUE ;
232     
233     for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
234       
235       if(phosNotFound){
236         sdigitsBranch=(TBranch *) branches->At(ibranch) ;
237         if( (strcmp("PHOS",sdigitsBranch->GetName())==0 ) &&
238             (fSDigitsTitle.CompareTo(sdigitsBranch->GetTitle()) == 0) )
239           phosNotFound = kFALSE ;
240       }
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 ;
246       }
247     }
248
249     if(!(sdigitizerNotFound && phosNotFound)){
250       cout << "AliPHOSSdigitizer error:" << endl ;
251       cout << "Can not overwrite existing branches: do not write" << endl ;
252       return ;
253     }
254     
255     //Make (if necessary) branches    
256     char * file =0;
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()) ;
260     }
261     
262     TDirectory *cwd = gDirectory;
263     
264     //First list of sdigits
265     Int_t bufferSize = 32000 ;    
266     sdigitsBranch = gAlice->TreeS()->Branch("PHOS",&fSDigits,bufferSize);
267     sdigitsBranch->SetTitle(fSDigitsTitle.Data());
268     if (file) {
269       sdigitsBranch->SetFile(file);
270       TIter next( sdigitsBranch->GetListOfBranches());
271       while ((sdigitsBranch=(TBranch*)next())) {
272         sdigitsBranch->SetFile(file);
273       }   
274       cwd->cd();
275     } 
276       
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());
283     if (file) {
284       sdigitizerBranch->SetFile(file);
285       TIter next( sdigitizerBranch->GetListOfBranches());
286       while ((sdigitizerBranch=(TBranch*)next())) {
287         sdigitizerBranch->SetFile(file);
288       }   
289       cwd->cd();
290       delete file;
291     }
292
293     gAlice->TreeS()->Fill() ;
294     gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
295     
296     if(strstr(option,"deb"))
297       PrintSDigits(option) ;
298     
299   }
300   
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 ;
306     cout << endl ;
307   }
308   
309   
310 }
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 ;
317 }
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 ;
326   
327 }
328 //__________________________________________________________________
329 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const{
330   if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
331     return kTRUE ;
332   else
333     return kFALSE ;
334 }
335 //__________________________________________________________________
336 void AliPHOSSDigitizer::PrintSDigits(Option_t * option){
337   //Prints list of digits produced at the current pass of AliPHOSDigitizer
338   
339   cout << "AliPHOSSDigitizer: " << endl ;
340   cout << "       Number of entries in SDigits list  " << fSDigits->GetEntriesFast() << endl ;
341   cout << endl ;
342   
343   if(strstr(option,"all")){// print all digits
344     
345     //loop over digits
346     AliPHOSDigit * digit;
347     cout << "SDigit Id " << " Amplitude " <<  " Index "  <<  " Nprim " << " Primaries list " <<  endl;    
348     Int_t index ;
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() <<"  ";
354       
355       Int_t iprimary;
356       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
357         cout << setw(5)  <<  digit->GetPrimary(iprimary+1) << "  ";
358       cout << endl;      
359     }
360     
361   }
362 }