]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSSDigitizer.cxx
comments added
[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 {
126   // Initialization can not be done in the default constructor
127
128   if(!fIsInitialized){
129
130     if(fHeadersFile.IsNull())
131       fHeadersFile="galice.root" ;
132
133     TFile * file = (TFile*) gROOT->GetFile(fHeadersFile.Data() ) ;
134     
135     //if file was not opened yet, read gAlice
136     if(file == 0){
137       file = new TFile(fHeadersFile.Data(),"update") ;
138       gAlice = (AliRun *) file->Get("gAlice") ;
139     }
140     
141     fHits    = new TClonesArray("AliPHOSHit",1000);
142     fSDigits = new TClonesArray("AliPHOSDigit",1000);
143     
144     // add Task to //root/Tasks folder
145     TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
146     roottasks->Add(this) ; 
147     
148     fIsInitialized = kTRUE ;
149   }
150 }
151 //____________________________________________________________________________
152 void AliPHOSSDigitizer::Exec(Option_t *option) 
153
154   // Collects all hits in the same active volume into digit
155   
156   if(!fIsInitialized)
157     Init() ;
158
159   if(strstr(option,"tim"))
160     gBenchmark->Start("PHOSSDigitizer");
161   
162   fNevents = (Int_t) gAlice->TreeE()->GetEntries() ; 
163   
164   Int_t ievent ;
165   for(ievent = 0; ievent < fNevents; ievent++){
166     gAlice->GetEvent(ievent) ;
167     gAlice->SetEvent(ievent) ;
168     
169     if(gAlice->TreeH()==0){
170       cout << "AliPHOSSDigitizer: There is no Hit Tree" << endl;
171       return ;
172     }
173     
174     //set address of the hits 
175     TBranch * branch = gAlice->TreeH()->GetBranch("PHOS");
176     if (branch) 
177       branch->SetAddress(&fHits);
178     else{
179       cout << "ERROR in AliPHOSSDigitizer: "<< endl ;
180       cout << "      no branch PHOS in TreeH"<< endl ;
181       cout << "      do nothing " << endl ;
182       return ;
183     }
184     
185     fSDigits->Clear();
186     Int_t nSdigits = 0 ;
187     
188     
189     //Now made SDigits from hits, for PHOS it is the same, so just copy    
190     Int_t itrack ;
191     for (itrack=0; itrack < gAlice->GetNtrack(); itrack++){
192       
193       //=========== Get the PHOS branch from Hits Tree for the Primary track itrack
194       branch->GetEntry(itrack,0);
195       
196       Int_t i;
197       for ( i = 0 ; i < fHits->GetEntries() ; i++ ) {
198         AliPHOSHit * hit = (AliPHOSHit*)fHits->At(i) ;
199
200         // Assign primary number only if contribution is significant
201         if( hit->GetEnergy() > fPrimThreshold)
202           new((*fSDigits)[nSdigits]) AliPHOSDigit( hit->GetPrimary(), hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
203         else
204           new((*fSDigits)[nSdigits]) AliPHOSDigit( -1               , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
205         
206         nSdigits++ ;  
207         
208       } 
209       
210     } // loop over tracks
211     
212     fSDigits->Sort() ;
213     
214     nSdigits = fSDigits->GetEntriesFast() ;
215     fSDigits->Expand(nSdigits) ;
216     
217     Int_t i ;
218     for (i = 0 ; i < nSdigits ; i++) { 
219       AliPHOSDigit * digit = (AliPHOSDigit *) fSDigits->At(i) ; 
220       digit->SetIndexInList(i) ;     
221     }
222
223     if(gAlice->TreeS() == 0)
224       gAlice->MakeTree("S") ;
225     
226     //check, if this branch already exits?
227     TBranch * sdigitsBranch = 0;
228     TBranch * sdigitizerBranch = 0;
229     
230     TObjArray * branches = gAlice->TreeS()->GetListOfBranches() ;
231     Int_t ibranch;
232     Bool_t phosNotFound = kTRUE ;
233     Bool_t sdigitizerNotFound = kTRUE ;
234     
235     for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
236       
237       if(phosNotFound){
238         sdigitsBranch=(TBranch *) branches->At(ibranch) ;
239         if( (strcmp("PHOS",sdigitsBranch->GetName())==0 ) &&
240             (fSDigitsTitle.CompareTo(sdigitsBranch->GetTitle()) == 0) )
241           phosNotFound = kFALSE ;
242       }
243       if(sdigitizerNotFound){
244         sdigitizerBranch = (TBranch *) branches->At(ibranch) ;
245         if( (strcmp(sdigitizerBranch->GetName(),"AliPHOSSDigitizer") == 0)&&
246             (fSDigitsTitle.CompareTo(sdigitizerBranch->GetTitle()) == 0) )
247           sdigitizerNotFound = kFALSE ;
248       }
249     }
250
251     if(!(sdigitizerNotFound && phosNotFound)){
252       cout << "AliPHOSSdigitizer error:" << endl ;
253       cout << "Can not overwrite existing branches: do not write" << endl ;
254       return ;
255     }
256     
257     //Make (if necessary) branches    
258     char * file =0;
259     if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
260       file = new char[strlen(gAlice->GetBaseFile())+20] ;
261       sprintf(file,"%s/PHOS.SDigits.root",gAlice->GetBaseFile()) ;
262     }
263     
264     TDirectory *cwd = gDirectory;
265     
266     //First list of sdigits
267     Int_t bufferSize = 32000 ;    
268     sdigitsBranch = gAlice->TreeS()->Branch("PHOS",&fSDigits,bufferSize);
269     sdigitsBranch->SetTitle(fSDigitsTitle.Data());
270     if (file) {
271       sdigitsBranch->SetFile(file);
272       TIter next( sdigitsBranch->GetListOfBranches());
273       TBranch * subbr;
274       while ((subbr=(TBranch*)next())) {
275         subbr->SetFile(file);
276       }   
277       cwd->cd();
278     } 
279       
280     //second - SDigitizer
281     Int_t splitlevel = 0 ;
282     AliPHOSSDigitizer * sd = this ;
283     sdigitizerBranch = gAlice->TreeS()->Branch("AliPHOSSDigitizer","AliPHOSSDigitizer",
284                                                &sd,bufferSize,splitlevel); 
285     sdigitizerBranch->SetTitle(fSDigitsTitle.Data());
286     if (file) {
287       sdigitizerBranch->SetFile(file);
288       TIter next( sdigitizerBranch->GetListOfBranches());
289       TBranch * subbr ;
290       while ((subbr=(TBranch*)next())) {
291         subbr->SetFile(file);
292       }   
293       cwd->cd();
294       delete file;
295     }
296
297     sdigitsBranch->Fill() ;
298     sdigitizerBranch->Fill() ;
299     gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
300     
301     if(strstr(option,"deb"))
302       PrintSDigits(option) ;
303     
304   }
305   
306   if(strstr(option,"tim")){
307     gBenchmark->Stop("PHOSSDigitizer");
308     cout << "AliPHOSSDigitizer:" << endl ;
309     cout << "   took " << gBenchmark->GetCpuTime("PHOSSDigitizer") << " seconds for SDigitizing " 
310          <<  gBenchmark->GetCpuTime("PHOSSDigitizer")/fNevents << " seconds per event " << endl ;
311     cout << endl ;
312   }
313   
314   
315 }
316 //__________________________________________________________________
317 void AliPHOSSDigitizer::SetSDigitsBranch(const char * title )
318 {
319   // Setting title to branch SDigits 
320   if(!fSDigitsTitle.IsNull())
321     cout << "AliPHOSSdigitizer: changing SDigits file from " <<fSDigitsTitle.Data() << " to " << title << endl ;
322   fSDigitsTitle=title ;
323 }
324 //__________________________________________________________________
325 void AliPHOSSDigitizer::Print(Option_t* option)const
326 {
327   // Prints parameters of SDigitizer
328   cout << "------------------- "<< GetName() << " -------------" << endl ;
329   cout << "   Writing SDigitis to branch with title  " << fSDigitsTitle.Data() << endl ;
330   cout << "   with digitization parameters  A = " << fA << endl ;
331   cout << "                                 B = " << fB << endl ;
332   cout << "   Threshold for Primary assignment= " << fPrimThreshold << endl ; 
333   cout << "---------------------------------------------------"<<endl ;
334   
335 }
336 //__________________________________________________________________
337 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
338 {
339   // SDigitizer are identical if the same threshold is in use
340   if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
341     return kTRUE ;
342   else
343     return kFALSE ;
344 }
345 //__________________________________________________________________
346 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
347 {
348   // Prints list of digits produced in the current pass of AliPHOSDigitizer
349   
350   cout << "AliPHOSSDigitizer: " << endl ;
351   cout << "       Number of entries in SDigits list  " << fSDigits->GetEntriesFast() << endl ;
352   cout << endl ;
353   
354   if(strstr(option,"all")){// print all digits
355     
356     //loop over digits
357     AliPHOSDigit * digit;
358     cout << "SDigit Id " << " Amplitude " <<  " Index "  <<  " Nprim " << " Primaries list " <<  endl;    
359     Int_t index ;
360     for (index = 0 ; index < fSDigits->GetEntries() ; index++) {
361       digit = (AliPHOSDigit * )  fSDigits->At(index) ;
362       cout << setw(8)  <<  digit->GetId() << " "  <<    setw(3)  <<  digit->GetAmp() <<   "  "  
363            << setw(6)  <<  digit->GetIndexInList() << "  "   
364            << setw(5)  <<  digit->GetNprimary() <<"  ";
365       
366       Int_t iprimary;
367       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
368         cout << setw(5)  <<  digit->GetPrimary(iprimary+1) << "  ";
369       cout << endl;      
370     }
371     
372   }
373 }