]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSSDigitizer.cxx
libHBTAnalysis renamed to libHBTAN
[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 "AliPHOSGeometry.h"
67 #include "AliPHOSGetter.h"
68 #include "AliPHOSHit.h"
69 #include "AliPHOSSDigitizer.h"
70
71
72 ClassImp(AliPHOSSDigitizer)
73
74            
75 //____________________________________________________________________________ 
76   AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("","") 
77 {
78   // ctor
79   fA             = 0;
80   fB             = 10000000.;
81   fPrimThreshold = 0.01 ;
82   fSDigitsInRun  = 0 ; 
83 }
84
85 //____________________________________________________________________________ 
86 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * headerFile, const char * sDigitsTitle):TTask(sDigitsTitle, headerFile)
87 {
88   // ctor
89   fA             = 0;
90   fB             = 10000000.;
91   fPrimThreshold = 0.01 ;
92   fSDigitsInRun  = 0 ; 
93   Init();
94 }
95
96
97 //____________________________________________________________________________ 
98 void AliPHOSSDigitizer::Init()
99 {
100   // Initialization: open root-file, allocate arrays for hits and sdigits,
101   // attach task SDigitizer to the list of PHOS tasks
102   // 
103   // Initialization can not be done in the default constructor
104   //============================================================= YS
105   //  The initialisation is now done by AliPHOSGetter
106   
107   if( strcmp(GetTitle(), "") == 0 )
108     SetTitle("galice.root") ;
109   
110   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName()) ;     
111   if ( gime == 0 ) {
112     cerr << "ERROR: AliPHOSSDigitizer::Init -> Could not obtain the Getter object !" << endl ; 
113     return ;
114   } 
115   
116   gime->PostSDigits( GetName(), GetTitle() ) ; 
117
118   TString sdname(GetName() );
119   sdname.Append(":") ;
120   sdname.Append(GetTitle() ) ;
121   SetName(sdname) ;
122   gime->PostSDigitizer(this) ;
123  
124   // create a folder on the white board //YSAlice/WhiteBoard/SDigits/PHOS/headerFile/sdigitsTitle
125     
126 }
127
128 //____________________________________________________________________________
129 void AliPHOSSDigitizer::Exec(Option_t *option) 
130
131   // Collects all hits in the same active volume into digit
132   if( strcmp(GetName(), "") == 0 )
133     Init() ;
134   
135   if (strstr(option, "print") ) {
136     Print("") ; 
137     return ; 
138   }
139
140   if(strstr(option,"tim"))
141     gBenchmark->Start("PHOSSDigitizer");
142
143   //Check, if this branch already exits
144   gAlice->GetEvent(0) ;
145   if(gAlice->TreeS() ) {
146     TObjArray * lob = static_cast<TObjArray*>(gAlice->TreeS()->GetListOfBranches()) ;
147     TIter next(lob) ; 
148     TBranch * branch = 0 ;  
149     Bool_t phosfound = kFALSE, sdigitizerfound = kFALSE ; 
150     
151     while ( (branch = (static_cast<TBranch*>(next()))) && (!phosfound || !sdigitizerfound) ) {
152       if ( (strcmp(branch->GetName(), "PHOS")==0) && (strcmp(branch->GetTitle(), GetName())==0) ) 
153         phosfound = kTRUE ;
154       
155       else if ( (strcmp(branch->GetName(), "AliPHOSSDigitizer")==0) && (strcmp(branch->GetTitle(), GetName())==0) ) 
156         sdigitizerfound = kTRUE ; 
157     }
158     
159     if ( phosfound || sdigitizerfound ) {
160       cerr << "WARNING: AliPHOSSDigitizer::Exec -> SDigits and/or SDigitizer branch with name " << GetName() 
161            << " already exits" << endl ;
162       return ; 
163     }   
164   }  
165
166   TString sdname(GetName()) ;
167   sdname.Remove(sdname.Index(GetTitle())-1) ;
168     
169   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
170
171   Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ; 
172   
173   Int_t ievent ;
174   for(ievent = 0; ievent < nevents; ievent++){
175     gime->Event(ievent,"H") ;
176     const TClonesArray * hits = gime->Hits() ;
177     TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
178     sdigits->Clear();
179     Int_t nSdigits = 0 ;
180     
181     
182     //Now make SDigits from hits, for PHOS it is the same, so just copy    
183
184 //******************** CHECK HERE
185 //YS DOES NOT UNDERSTAND THE NEED FOR THE FOLLOWING LOOP
186 //     Int_t itrack ;
187 //     for (itrack=0; itrack < gAlice->GetNtrack(); itrack++){
188 //       //=========== Get the PHOS branch from Hits Tree for the Primary track itrack
189 //       gime->Track(itrack) ;
190       Int_t i;
191       for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
192         AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
193         // Assign primary number only if contribution is significant
194         
195         if( hit->GetEnergy() > fPrimThreshold)
196           new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
197                                                   Digitize(hit->GetEnergy()), hit->GetTime()) ;
198         else
199           new((*sdigits)[nSdigits]) AliPHOSDigit( -1              , hit->GetId(), 
200                                                    Digitize(hit->GetEnergy()), hit->GetTime()) ;
201         nSdigits++ ;    
202         
203       }
204 //    } // loop over tracks
205     
206     sdigits->Sort() ;
207     
208     nSdigits = sdigits->GetEntriesFast() ;
209     fSDigitsInRun += nSdigits ;  
210     sdigits->Expand(nSdigits) ;
211 //    Int_t i ;
212     for (i = 0 ; i < nSdigits ; i++) { 
213       AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ; 
214       digit->SetIndexInList(i) ;     
215     }
216
217     if(gAlice->TreeS() == 0)
218       gAlice->MakeTree("S") ;
219     
220     //Make (if necessary) branches    
221     char * file =0;
222     if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
223       file = new char[strlen(gAlice->GetBaseFile())+20] ;
224       sprintf(file,"%s/PHOS.SDigits.root",gAlice->GetBaseFile()) ;
225     }
226     
227     TDirectory *cwd = gDirectory;
228     
229     //First list of sdigits
230     Int_t bufferSize = 32000 ;    
231     TBranch * sdigitsBranch = gAlice->TreeS()->Branch("PHOS",&sdigits,bufferSize);
232     sdigitsBranch->SetTitle(sdname);
233     if (file) {
234       sdigitsBranch->SetFile(file);
235       TIter next( sdigitsBranch->GetListOfBranches());
236       TBranch * subbr;
237       while ((subbr=static_cast<TBranch*>(next()))) {
238         subbr->SetFile(file);
239       }   
240       cwd->cd();
241     } 
242       
243     //Next - SDigitizer
244     Int_t splitlevel = 0 ;
245     AliPHOSSDigitizer * sd = this ;
246     TBranch * sdigitizerBranch = gAlice->TreeS()->Branch("AliPHOSSDigitizer","AliPHOSSDigitizer",
247                                                &sd,bufferSize,splitlevel); 
248     sdigitizerBranch->SetTitle(sdname);
249     if (file) {
250       sdigitizerBranch->SetFile(file);
251       TIter next( sdigitizerBranch->GetListOfBranches());
252       TBranch * subbr ;
253       while ((subbr=static_cast<TBranch*>(next()))) {
254         subbr->SetFile(file);
255       }   
256       cwd->cd();
257       delete file;
258     }
259
260     sdigitsBranch->Fill() ;
261     sdigitizerBranch->Fill() ;
262     gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
263     
264     if(strstr(option,"deb"))
265       PrintSDigits(option) ;
266     
267   }
268   
269   if(strstr(option,"tim")){
270     gBenchmark->Stop("PHOSSDigitizer");
271     cout << "AliPHOSSDigitizer:" << endl ;
272     cout << "   took " << gBenchmark->GetCpuTime("PHOSSDigitizer") << " seconds for SDigitizing " 
273          <<  gBenchmark->GetCpuTime("PHOSSDigitizer")/nevents << " seconds per event " << endl ;
274     cout << endl ;
275   }
276   
277   
278 }
279 //__________________________________________________________________
280 void AliPHOSSDigitizer::SetSDigitsBranch(const char * title )
281 {
282   // Setting title to branch SDigits 
283
284   TString stitle(title) ;
285
286   // check if branch with title already exists
287   TBranch * sdigitsBranch    = 
288     static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("PHOS")) ; 
289   TBranch * sdigitizerBranch =  
290     static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("AliPHOSSDigitizer")) ;
291   const char * sdigitsTitle    = sdigitsBranch ->GetTitle() ;  
292   const char * sdigitizerTitle = sdigitizerBranch ->GetTitle() ;
293   if ( stitle.CompareTo(sdigitsTitle)==0 || stitle.CompareTo(sdigitizerTitle)==0 ){
294     cerr << "ERROR: AliPHOSSdigitizer::SetSDigitsBranch -> Cannot overwrite existing branch with title " << title << endl ;
295     return ;
296   }
297   
298   cout << "AliPHOSSdigitizer::SetSDigitsBranch -> Changing SDigits file from " << GetName() << " to " << title << endl ;
299
300   SetName(title) ; 
301     
302   // Post to the WhiteBoard
303   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
304   gime->PostSDigits( title, GetTitle()) ; 
305 }
306
307 //__________________________________________________________________
308 void AliPHOSSDigitizer::Print(Option_t* option)const
309 {
310   // Prints parameters of SDigitizer
311   cout << "------------------- "<< GetName() << " -------------" << endl ;
312   cout << "   Writing SDigits to branch with title  " << GetName() << endl ;
313   cout << "   with digitization parameters  A = " << fA << endl ;
314   cout << "                                 B = " << fB << endl ;
315   cout << "   Threshold for Primary assignment= " << fPrimThreshold << endl ; 
316   cout << "---------------------------------------------------"<<endl ;
317   
318 }
319 //__________________________________________________________________
320 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
321 {
322   // Equal operator.
323   // SDititizers are equal if their pedestal, slope and threshold are equal
324
325   if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
326     return kTRUE ;
327   else
328     return kFALSE ;
329 }
330 //__________________________________________________________________
331 //__________________________________________________________________
332 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
333 {
334   // Prints list of digits produced in the current pass of AliPHOSDigitizer
335
336
337   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
338   TString sdname(GetName()) ;
339   sdname.Remove(sdname.Index(GetTitle())-1) ;
340   const TClonesArray * sdigits = gime->SDigits(sdname.Data()) ; 
341
342   cout << "AliPHOSSDigitiser: event " << gAlice->GetEvNumber() << endl ;
343   cout << "      Number of entries in SDigits list " << sdigits->GetEntriesFast() << endl ;
344   cout << endl ;
345   if(strstr(option,"all")||strstr(option,"EMC")){
346     
347     //loop over digits
348     AliPHOSDigit * digit;
349     cout << "EMC sdigits " << endl ;
350     cout << "Digit Id    Amplitude     Index     Nprim  Primaries list " <<  endl;      
351     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
352     Int_t index ;
353     for (index = 0 ; (index < sdigits->GetEntriesFast()) && 
354          (((AliPHOSDigit * )  sdigits->At(index))->GetId() <= maxEmc) ; index++) {
355       digit = (AliPHOSDigit * )  sdigits->At(index) ;
356       if(digit->GetNprimary() == 0) continue;
357       cout << setw(6)  <<  digit->GetId() << "   "  <<  setw(10)  <<  digit->GetAmp() <<   "    "  
358            << setw(6)  <<  digit->GetIndexInList() << "    "   
359            << setw(5)  <<  digit->GetNprimary() <<"    ";
360       
361       Int_t iprimary;
362       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
363         cout << setw(5)  <<  digit->GetPrimary(iprimary+1) << "    ";
364       cout << endl;      
365     }    
366     cout << endl;
367   }
368
369   if(strstr(option,"all")||strstr(option,"CPV")){
370     
371     //loop over CPV digits
372     AliPHOSDigit * digit;
373     cout << "CPV sdigits " << endl ;
374     cout << "Digit Id    Amplitude     Index     Nprim  Primaries list " <<  endl;      
375     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
376     Int_t index ;
377     for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
378       digit = (AliPHOSDigit * )  sdigits->At(index) ;
379       if(digit->GetId() > maxEmc){
380         cout << setw(6)  <<  digit->GetId() << "   "  <<        setw(10)  <<  digit->GetAmp() <<   "    "  
381              << setw(6)  <<  digit->GetIndexInList() << "    "   
382              << setw(5)  <<  digit->GetNprimary() <<"    ";
383         
384         Int_t iprimary;
385         for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
386           cout << setw(5)  <<  digit->GetPrimary(iprimary+1) << "    ";
387         cout << endl;    
388       }    
389     }
390   }
391
392 }
393
394 //____________________________________________________________________________ 
395 void AliPHOSSDigitizer::UseHitsFrom(const char * filename)
396 {
397   SetTitle(filename) ; 
398   Init() ; 
399 }