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