]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSSDigitizer.cxx
Mostly minor style modifications to be ready for cloning with EMCAL
[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
17 /* $Id$ */
18
19 //_________________________________________________________________________
20 // This is a TTask that makes SDigits out of Hits
21 // The name of the TTask is also the title of the branch that will contain 
22 // the created SDigits
23 // The title of the TTAsk is the name of the file that contains the hits from
24 // which the SDigits are created
25 // A Summable Digits is the sum of all hits originating 
26 // from one primary in one active cell
27 // A threshold for assignment of the primary to SDigit is applied 
28 // SDigits are written to TreeS, branch "PHOS"
29 // AliPHOSSDigitizer with all current parameters is written 
30 // to TreeS branch "AliPHOSSDigitizer".
31 // Both branches have the same title. If necessary one can produce 
32 // another set of SDigits with different parameters. Two versions
33 // can be distunguished using titles of the branches.
34 // User case:
35 //  root [0] AliPHOSSDigitizer * s = new AliPHOSSDigitizer("galice.root")
36 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
37 //  root [1] s->ExecuteTask()
38 //             // Makes SDigitis for all events stored in galice.root
39 //  root [2] s->SetPedestalParameter(0.001)
40 //             // One can change parameters of digitization
41 // root [3] s->SetSDigitsBranch("Pedestal 0.001")
42 //             // and write them into the new branch
43 // root [4] s->ExecuteTask("deb all tim")
44 //             // available parameters:
45 //             deb - print # of produced SDigitis
46 //             deb all  - print # and list of produced SDigits
47 //             tim - print benchmarking information
48 //
49 //*-- Author :  Dmitri Peressounko (SUBATECH & KI) 
50 //////////////////////////////////////////////////////////////////////////////
51
52
53 // --- ROOT system ---
54 #include "TFile.h"
55 #include "TTask.h"
56 #include "TTree.h"
57 #include "TSystem.h"
58 #include "TROOT.h"
59 #include "TFolder.h"
60 #include "TBenchmark.h"
61 #include "TGeometry.h"
62
63 // --- Standard library ---
64 #include <iomanip.h>
65
66 // --- AliRoot header files ---
67 #include "AliRun.h"
68 #include "AliHeader.h"
69 #include "AliPHOSDigit.h"
70 #include "AliPHOSGeometry.h"
71 #include "AliPHOSGetter.h"
72 #include "AliPHOSHit.h"
73 #include "AliPHOSSDigitizer.h"
74
75 ClassImp(AliPHOSSDigitizer)
76
77            
78 //____________________________________________________________________________ 
79   AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("","") 
80 {
81   // ctor
82   InitParameters() ;
83   fDefaultInit = kTRUE ; 
84 }
85
86 //____________________________________________________________________________ 
87 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * headerFile, const char * sDigitsTitle, const Bool_t toSplit):
88 TTask(sDigitsTitle, headerFile)
89 {
90   // ctor
91   InitParameters() ; 
92   fToSplit = toSplit ;
93   Init();
94   fDefaultInit = kFALSE ; 
95 }
96
97 //____________________________________________________________________________ 
98 AliPHOSSDigitizer::~AliPHOSSDigitizer()
99 {
100   // dtor
101   
102   fSplitFile = 0 ; 
103 }
104
105 //____________________________________________________________________________ 
106 void AliPHOSSDigitizer::Init()
107 {
108   // Initialization: open root-file, allocate arrays for hits and sdigits,
109   // attach task SDigitizer to the list of PHOS tasks
110   // 
111   // Initialization can not be done in the default constructor
112   //============================================================= YS
113   //  The initialisation is now done by AliPHOSGetter
114   
115   if( strcmp(GetTitle(), "") == 0 )
116     SetTitle("galice.root") ;
117
118   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName(),fToSplit) ;  
119   if ( gime == 0 ) {
120     cerr << "ERROR: AliPHOSSDigitizer::Init -> Could not obtain the Getter object !" 
121          << endl ; 
122     return ;
123   } 
124   
125   gime->PostSDigits( GetName(), GetTitle() ) ; 
126
127   fSplitFile = 0 ;
128   if(fToSplit){
129     // construct the name of the file as /path/PHOS.SDigits.root
130     // First - extract full path if necessary
131     TString sDigitsFileName(GetTitle()) ;
132     Ssiz_t islash = sDigitsFileName.Last('/') ;
133     if(islash<sDigitsFileName.Length())
134       sDigitsFileName.Remove(islash+1,sDigitsFileName.Length()) ;
135     else
136       sDigitsFileName="" ;
137     // Next - append the file name 
138     sDigitsFileName+="PHOS.SDigits." ;
139     if((strcmp(GetName(),"Default")!=0)&&(strcmp(GetName(),"")!=0)){
140       sDigitsFileName+=GetName() ;
141       sDigitsFileName+="." ;
142     }
143     sDigitsFileName+="root" ;
144     // Finally - check if the file already opened or open the file
145     fSplitFile = static_cast<TFile*>(gROOT->GetFile(sDigitsFileName.Data()));   
146     if(!fSplitFile)
147       fSplitFile =  TFile::Open(sDigitsFileName.Data(),"update") ;
148   }
149
150   TString sdname(GetName() );
151   sdname.Append(":") ;
152   sdname.Append(GetTitle() ) ;
153   SetName(sdname) ;
154   gime->PostSDigitizer(this) ;
155 }
156
157 //____________________________________________________________________________ 
158 void AliPHOSSDigitizer::InitParameters()
159
160   fA             = 0;
161   fB             = 10000000.;
162   fPrimThreshold = 0.01 ;
163   fSDigitsInRun  = 0 ;
164   fSplitFile     = 0 ; 
165   fToSplit       = kFALSE ;
166 }
167
168 //____________________________________________________________________________
169 void AliPHOSSDigitizer::Exec(Option_t *option) 
170
171   // Collects all hits in the same active volume into digit
172
173   if( strcmp(GetName(), "") == 0 )
174     Init() ;
175   
176   if (strstr(option, "print") ) {
177     Print("") ; 
178     return ; 
179   }
180
181   if(strstr(option,"tim"))
182     gBenchmark->Start("PHOSSDigitizer");
183   
184   //Check, if this branch already exits
185   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
186   if(gime->BranchExists("SDigits") ) 
187     return;   
188
189   TString sdname(GetName()) ;
190   sdname.Remove(sdname.Index(GetTitle())-1) ;
191    
192   Int_t nevents = gime->MaxEvent() ; 
193   Int_t ievent ;
194   for(ievent = 0; ievent < nevents; ievent++){
195     gime->Event(ievent,"H") ;
196     const TClonesArray * hits = gime->Hits() ;
197     TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
198     sdigits->Clear();
199     Int_t nSdigits = 0 ;
200     
201     //Now make SDigits from hits, for PHOS it is the same, so just copy    
202     Int_t nPrim =  static_cast<Int_t>((gAlice->TreeH())->GetEntries()) ; 
203     // Attention nPrim is the number of primaries tracked by Geant 
204     // and this number could be different to the number of Primaries in TreeK;
205     Int_t iprim ;
206     for (iprim = 0 ; iprim < nPrim ; iprim ++) { 
207       //=========== Get the PHOS branch from Hits Tree for the Primary iprim
208       gime->Track(iprim) ;
209       Int_t i;
210       for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
211         AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
212         // Assign primary number only if contribution is significant
213         
214         if( hit->GetEnergy() > fPrimThreshold)
215           new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
216                                                  Digitize(hit->GetEnergy()), hit->GetTime()) ;
217         else
218           new((*sdigits)[nSdigits]) AliPHOSDigit( -1              , hit->GetId(), 
219                                                   Digitize(hit->GetEnergy()), hit->GetTime()) ;
220         nSdigits++ ;    
221         
222       }
223       
224     } // loop over iprim
225     
226     sdigits->Sort() ;
227     
228     nSdigits = sdigits->GetEntriesFast() ;
229     fSDigitsInRun += nSdigits ;  
230     sdigits->Expand(nSdigits) ;
231     
232     Int_t i ;
233     for (i = 0 ; i < nSdigits ; i++) { 
234       AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ; 
235       digit->SetIndexInList(i) ;     
236     }
237
238     //Now write SDigits
239     
240     if((gAlice->TreeS() == 0)|| (fSplitFile))  
241       gAlice->MakeTree("S", fSplitFile);
242     
243     if(fSplitFile)
244       fSplitFile->cd() ;
245
246     //First list of sdigits
247     Int_t bufferSize = 32000 ;
248     TBranch * sdigitsBranch = gAlice->TreeS()->Branch("PHOS",&sdigits,bufferSize);
249     sdigitsBranch->SetTitle(sdname);
250     
251     //Next - SDigitizer
252     Int_t splitlevel = 0 ;
253     AliPHOSSDigitizer * sd = this ;
254     TBranch * sdigitizerBranch = gAlice->TreeS()->Branch("AliPHOSSDigitizer","AliPHOSSDigitizer",
255                                                          &sd,bufferSize,splitlevel);
256     sdigitizerBranch->SetTitle(sdname);
257     
258     sdigitsBranch->Fill() ;
259     sdigitizerBranch->Fill() ;
260
261     gAlice->TreeS()->AutoSave() ;
262         
263     if(strstr(option,"deb"))
264       PrintSDigits(option) ;
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 void AliPHOSSDigitizer::SetSDigitsBranch(const char * title )
278 {
279   // Setting title to branch SDigits 
280
281   TString stitle(title) ;
282
283   // check if branch with title already exists
284   TBranch * sdigitsBranch    = 
285     static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("PHOS")) ; 
286   TBranch * sdigitizerBranch =  
287     static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("AliPHOSSDigitizer")) ;
288   const char * sdigitsTitle    = sdigitsBranch ->GetTitle() ;  
289   const char * sdigitizerTitle = sdigitizerBranch ->GetTitle() ;
290   if ( stitle.CompareTo(sdigitsTitle)==0 || stitle.CompareTo(sdigitizerTitle)==0 ){
291     cerr << "ERROR: AliPHOSSdigitizer::SetSDigitsBranch -> Cannot overwrite existing branch with title " << title << endl ;
292     return ;
293   }
294   
295   cout << "AliPHOSSdigitizer::SetSDigitsBranch -> Changing SDigits file from " << GetName() << " to " << title << endl ;
296
297   SetName(title) ; 
298     
299   // Post to the WhiteBoard
300   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
301   gime->PostSDigits( title, GetTitle()) ; 
302 }
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 //__________________________________________________________________
319 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
320 {
321   // Equal operator.
322   // SDititizers are equal if their pedestal, slope and threshold are equal
323
324   if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
325     return kTRUE ;
326   else
327     return kFALSE ;
328 }
329
330 //__________________________________________________________________
331 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
332 {
333   // Prints list of digits produced in the current pass of AliPHOSDigitizer
334
335
336   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
337   TString sdname(GetName()) ;
338   sdname.Remove(sdname.Index(GetTitle())-1) ;
339   const TClonesArray * sdigits = gime->SDigits(sdname.Data()) ; 
340
341   cout << "AliPHOSSDigitiser: event " << gAlice->GetEvNumber() << endl ;
342   cout << "      Number of entries in SDigits list " << sdigits->GetEntriesFast() << endl ;
343   cout << endl ;
344   if(strstr(option,"all")||strstr(option,"EMC")){
345     
346     //loop over digits
347     AliPHOSDigit * digit;
348     cout << "EMC sdigits " << endl ;
349     cout << "Digit Id    Amplitude     Index     Nprim  Primaries list " <<  endl;      
350     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
351     Int_t index ;
352     for (index = 0 ; (index < sdigits->GetEntriesFast()) && 
353          ((dynamic_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
354       digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
355       if(digit->GetNprimary() == 0) continue;
356       cout << setw(6)  <<  digit->GetId() << "   "  <<  setw(10)  <<  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     cout << endl;
366   }
367
368   if(strstr(option,"all")||strstr(option,"CPV")){
369     
370     //loop over CPV digits
371     AliPHOSDigit * digit;
372     cout << "CPV sdigits " << endl ;
373     cout << "Digit Id    Amplitude     Index     Nprim  Primaries list " <<  endl;      
374     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
375     Int_t index ;
376     for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
377       digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
378       if(digit->GetId() > maxEmc){
379         cout << setw(6)  <<  digit->GetId() << "   "  <<        setw(10)  <<  digit->GetAmp() <<   "    "  
380              << setw(6)  <<  digit->GetIndexInList() << "    "   
381              << setw(5)  <<  digit->GetNprimary() <<"    ";
382         
383         Int_t iprimary;
384         for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
385           cout << setw(5)  <<  digit->GetPrimary(iprimary+1) << "    ";
386         cout << endl;    
387       }    
388     }
389   }
390
391 }
392
393 //____________________________________________________________________________ 
394 void AliPHOSSDigitizer::UseHitsFrom(const char * filename)
395 {
396   SetTitle(filename) ; 
397   Init() ; 
398 }