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