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