]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSSDigitizer.cxx
A lot of mods to put PHOS in the official TFolder/TTask structure
[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 // The name of the TTask is also the title of the branch that will contain 
21 // the created SDigits
22 // The title of the TTAsk is the name of the file that contains the hits from
23 // which the SDigits are created
24 // A Summable Digits is the sum of all hits originating 
25 // from one primary in one active cell
26 // A threshold for assignment of the primary to SDigit is applied 
27 // SDigits are written to TreeS, branch "PHOS"
28 // AliPHOSSDigitizer with all current parameters is written 
29 // to TreeS branch "AliPHOSSDigitizer".
30 // Both branches have the same title. If necessary one can produce 
31 // another set of SDigits with different parameters. Two versions
32 // can be distunguished using titles of the branches.
33 // User case:
34 //  root [0] AliPHOSSDigitizer * s = new AliPHOSSDigitizer("galice.root")
35 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
36 //  root [1] s->ExecuteTask()
37 //             // Makes SDigitis for all events stored in galice.root
38 //  root [2] s->SetPedestalParameter(0.001)
39 //             // One can change parameters of digitization
40 // root [3] s->SetSDigitsBranch("Pedestal 0.001")
41 //             // and write them into the new branch
42 // root [4] s->ExecuteTask("deb all tim")
43 //             // available parameters:
44 //             deb - print # of produced SDigitis
45 //             deb all  - print # and list of produced SDigits
46 //             tim - print benchmarking information
47 //
48 //*-- Author :  Dmitri Peressounko (SUBATECH & KI) 
49 //////////////////////////////////////////////////////////////////////////////
50
51
52 // --- ROOT system ---
53 #include "TFile.h"
54 #include "TTask.h"
55 #include "TTree.h"
56 #include "TSystem.h"
57 #include "TROOT.h"
58 #include "TFolder.h"
59 #include "TBenchmark.h"
60 // --- Standard library ---
61 #include <iomanip.h>
62
63 // --- AliRoot header files ---
64 #include "AliRun.h"
65 #include "AliPHOSDigit.h"
66 #include "AliPHOSGetter.h"
67 #include "AliPHOSHit.h"
68 #include "AliPHOSSDigitizer.h"
69
70
71 ClassImp(AliPHOSSDigitizer)
72
73            
74 //____________________________________________________________________________ 
75   AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("","") 
76 {
77   // ctor
78   fA             = 0;
79   fB             = 10000000.;
80   fPrimThreshold = 0.01 ;
81   fSDigitsInRun  = 0 ; 
82 }
83
84 //____________________________________________________________________________ 
85 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * headerFile, const char * sDigitsTitle):TTask(sDigitsTitle, headerFile)
86 {
87   // ctor
88   fA             = 0;
89   fB             = 10000000.;
90   fPrimThreshold = 0.01 ;
91   fSDigitsInRun  = 0 ; 
92   Init();
93 }
94
95
96 //____________________________________________________________________________ 
97 void AliPHOSSDigitizer::Init()
98 {
99   // Initialization: open root-file, allocate arrays for hits and sdigits,
100   // attach task SDigitizer to the list of PHOS tasks
101   // 
102   // Initialization can not be done in the default constructor
103   //============================================================= YS
104   //  The initialisation is now done by AliPHOSGetter
105   
106   if( strcmp(GetTitle(), "") == 0 )
107     SetTitle("galice.root") ;
108   
109   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName()) ;     
110   if ( gime == 0 ) {
111     cerr << "ERROR: AliPHOSSDigitizer::Init -> Could not obtain the Getter object !" << endl ; 
112     return ;
113   } 
114   
115   gime->PostSDigits( GetName(), GetTitle() ) ; 
116
117   TString sdname(GetName() );
118   sdname.Append(":") ;
119   sdname.Append(GetTitle() ) ;
120   SetName(sdname) ;
121   gime->PostSDigitizer(this) ;
122  
123   // create a folder on the white board //YSAlice/WhiteBoard/SDigits/PHOS/headerFile/sdigitsTitle
124     
125 }
126
127 //____________________________________________________________________________
128 void AliPHOSSDigitizer::Exec(Option_t *option) 
129
130   // Collects all hits in the same active volume into digit
131   if( strcmp(GetName(), "") == 0 )
132     Init() ;
133   
134   if (strstr(option, "print") ) {
135     Print("") ; 
136     return ; 
137   }
138
139   if(strstr(option,"tim"))
140     gBenchmark->Start("PHOSSDigitizer");
141
142   //Check, if this branch already exits
143   gAlice->GetEvent(0) ;
144   if(gAlice->TreeS() ) {
145     TObjArray * lob = static_cast<TObjArray*>(gAlice->TreeS()->GetListOfBranches()) ;
146     TIter next(lob) ; 
147     TBranch * branch = 0 ;  
148     Bool_t phosfound = kFALSE, sdigitizerfound = kFALSE ; 
149     
150     while ( (branch = (static_cast<TBranch*>(next()))) && (!phosfound || !sdigitizerfound) ) {
151       if ( (strcmp(branch->GetName(), "PHOS")==0) && (strcmp(branch->GetTitle(), GetName())==0) ) 
152         phosfound = kTRUE ;
153       
154       else if ( (strcmp(branch->GetName(), "AliPHOSSDigitizer")==0) && (strcmp(branch->GetTitle(), GetName())==0) ) 
155         sdigitizerfound = kTRUE ; 
156     }
157     
158     if ( phosfound || sdigitizerfound ) {
159       cerr << "WARNING: AliPHOSSDigitizer::Exec -> SDigits and/or SDigitizer branch with name " << GetName() 
160            << " already exits" << endl ;
161       return ; 
162     }   
163   }  
164
165   TString sdname(GetName()) ;
166   sdname.Remove(sdname.Index(GetTitle())-1) ;
167     
168   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
169
170   Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ; 
171   
172   Int_t ievent ;
173   for(ievent = 0; ievent < nevents; ievent++){
174     gime->Event(ievent,"H") ;
175
176     TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
177     sdigits->Clear();
178     Int_t nSdigits = 0 ;
179     
180     
181     //Now make 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       gime->Track(itrack) ;
187       TClonesArray * hits = gime->Hits() ;
188       Int_t i;
189       for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
190         AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
191         // Assign primary number only if contribution is significant
192         
193         if( hit->GetEnergy() > fPrimThreshold)
194           new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
195                                                   Digitize(hit->GetEnergy()), hit->GetTime()) ;
196         else
197           new((*sdigits)[nSdigits]) AliPHOSDigit( -1              , hit->GetId(), 
198                                                    Digitize(hit->GetEnergy()), hit->GetTime()) ;
199         nSdigits++ ;    
200         
201       }
202     } // loop over tracks
203     
204     sdigits->Sort() ;
205     
206     nSdigits = sdigits->GetEntriesFast() ;
207     
208     sdigits->Expand(nSdigits) ;
209     Int_t i ;
210     for (i = 0 ; i < nSdigits ; i++) { 
211       AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ; 
212       digit->SetIndexInList(i) ;     
213     }
214
215     if(gAlice->TreeS() == 0)
216       gAlice->MakeTree("S") ;
217     
218     //Make (if necessary) branches    
219     char * file =0;
220     if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
221       file = new char[strlen(gAlice->GetBaseFile())+20] ;
222       sprintf(file,"%s/PHOS.SDigits.root",gAlice->GetBaseFile()) ;
223     }
224     
225     TDirectory *cwd = gDirectory;
226     
227     //First list of sdigits
228     Int_t bufferSize = 32000 ;    
229     TBranch * sdigitsBranch = gAlice->TreeS()->Branch("PHOS",&sdigits,bufferSize);
230     sdigitsBranch->SetTitle(sdname);
231     if (file) {
232       sdigitsBranch->SetFile(file);
233       TIter next( sdigitsBranch->GetListOfBranches());
234       TBranch * subbr;
235       while ((subbr=static_cast<TBranch*>(next()))) {
236         subbr->SetFile(file);
237       }   
238       cwd->cd();
239     } 
240       
241     //Next - SDigitizer
242     Int_t splitlevel = 0 ;
243     AliPHOSSDigitizer * sd = this ;
244     TBranch * sdigitizerBranch = gAlice->TreeS()->Branch("AliPHOSSDigitizer","AliPHOSSDigitizer",
245                                                &sd,bufferSize,splitlevel); 
246     sdigitizerBranch->SetTitle(sdname);
247     if (file) {
248       sdigitizerBranch->SetFile(file);
249       TIter next( sdigitizerBranch->GetListOfBranches());
250       TBranch * subbr ;
251       while ((subbr=static_cast<TBranch*>(next()))) {
252         subbr->SetFile(file);
253       }   
254       cwd->cd();
255       delete file;
256     }
257
258     sdigitsBranch->Fill() ;
259     sdigitizerBranch->Fill() ;
260     gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
261     
262     if(strstr(option,"deb"))
263       PrintSDigits(option) ;
264     
265   }
266   
267   if(strstr(option,"tim")){
268     gBenchmark->Stop("PHOSSDigitizer");
269     cout << "AliPHOSSDigitizer:" << endl ;
270     cout << "   took " << gBenchmark->GetCpuTime("PHOSSDigitizer") << " seconds for SDigitizing " 
271          <<  gBenchmark->GetCpuTime("PHOSSDigitizer")/nevents << " seconds per event " << endl ;
272     cout << endl ;
273   }
274   
275   
276 }
277 //__________________________________________________________________
278 void AliPHOSSDigitizer::SetSDigitsBranch(const char * title )
279 {
280   // Setting title to branch SDigits 
281
282   TString stitle(title) ;
283
284   // check if branch with title already exists
285   TBranch * sdigitsBranch    = 
286     static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("PHOS")) ; 
287   TBranch * sdigitizerBranch =  
288     static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("AliPHOSSDigitizer")) ;
289   const char * sdigitsTitle    = sdigitsBranch ->GetTitle() ;  
290   const char * sdigitizerTitle = sdigitizerBranch ->GetTitle() ;
291   if ( stitle.CompareTo(sdigitsTitle)==0 || stitle.CompareTo(sdigitizerTitle)==0 ){
292     cerr << "ERROR: AliPHOSSdigitizer::SetSDigitsBranch -> Cannot overwrite existing branch with title " << title << endl ;
293     return ;
294   }
295   
296   cout << "AliPHOSSdigitizer::SetSDigitsBranch -> Changing SDigits file from " << GetName() << " to " << title << endl ;
297
298   SetName(title) ; 
299     
300   // Post to the WhiteBoard
301   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
302   gime->PostSDigits( title, GetTitle()) ; 
303 }
304
305 //__________________________________________________________________
306 void AliPHOSSDigitizer::Print(Option_t* option)const
307 {
308   // Prints parameters of SDigitizer
309   cout << "------------------- "<< GetName() << " -------------" << endl ;
310   cout << "   Writing SDigits to branch with title  " << GetName() << endl ;
311   cout << "   with digitization parameters  A = " << fA << endl ;
312   cout << "                                 B = " << fB << endl ;
313   cout << "   Threshold for Primary assignment= " << fPrimThreshold << endl ; 
314   cout << "---------------------------------------------------"<<endl ;
315   
316 }
317 //__________________________________________________________________
318 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
319 {
320   // Equal operator.
321   // SDititizers are equal if their pedestal, slope and threshold are equal
322
323   if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
324     return kTRUE ;
325   else
326     return kFALSE ;
327 }
328 //__________________________________________________________________
329 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
330 {
331   // Prints list of digits produced in the current pass of AliPHOSDigitizer
332
333   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
334   TString sdname(GetName()) ;
335   sdname.Remove(sdname.Index(GetTitle())-1) ;
336   TClonesArray * sdigits = gime->SDigits(sdname.Data()) ; 
337
338   cout << "AliPHOSSDigitizer: event "<<gAlice->GetEvNumber() << endl ;
339   cout << "       Number of entries in SDigits list  " << sdigits->GetEntriesFast() << endl ;
340   cout << endl ;
341
342   fSDigitsInRun +=  sdigits->GetEntriesFast() ; 
343   
344   if(strstr(option,"all")){// print all digits
345     
346     //loop over digits
347     AliPHOSDigit * digit;
348     cout << "SDigit Id " << " Amplitude " <<  " Index "  <<  " Nprim " << " Primaries list " <<  endl;    
349     Int_t index ;
350     for (index = 0 ; index < sdigits->GetEntries() ; index++) {
351       digit = (AliPHOSDigit * )  sdigits->At(index) ;
352       cout << setw(8)  <<  digit->GetId() << " "  <<    setw(3)  <<  digit->GetAmp() <<   "  "  
353            << setw(6)  <<  digit->GetIndexInList() << "  "   
354            << setw(5)  <<  digit->GetNprimary() <<"  ";
355       
356       Int_t iprimary;
357       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
358         cout << setw(5)  <<  digit->GetPrimary(iprimary+1) << "  ";
359       cout << endl;      
360     }
361     
362   }
363 }
364
365 //____________________________________________________________________________ 
366 void AliPHOSSDigitizer::UseHitsFrom(const char * filename)
367 {
368   SetTitle(filename) ; 
369   Init() ; 
370 }