]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSSDigitizer.cxx
More exact rounding function, but also much slower.
[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
76 ClassImp(AliPHOSSDigitizer)
77
78            
79 //____________________________________________________________________________ 
80   AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("","") {
81   // ctor
82   InitParameters() ;
83   fDefaultInit = kTRUE ; 
84 }
85
86 //____________________________________________________________________________ 
87 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * headerFile, const char * sDigitsTitle):TTask(sDigitsTitle, headerFile)
88 {
89   // ctor
90   InitParameters() ; 
91   Init();
92   fDefaultInit = kFALSE ; 
93 }
94
95 //____________________________________________________________________________ 
96 AliPHOSSDigitizer::~AliPHOSSDigitizer()
97 {
98   // dtor
99   // fDefaultInit = kTRUE if SDigitizer created by default ctor (to get just the parameters)
100
101   
102   if (!fDefaultInit) {
103     AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
104     
105     // remove the task from the folder list
106     gime->RemoveTask("S",GetName()) ;
107     
108     TString name(GetName()) ; 
109     if (! name.IsNull() ) 
110       name.Remove(name.Index(":")) ; 
111     
112     // remove the Hits from the folder list
113     gime->RemoveObjects("H",name) ;
114     
115     // remove the SDigits from the folder list
116     gime->RemoveObjects("S", name) ;
117     
118     // Delete gAlice
119     gime->CloseFile() ; 
120     
121   }
122   fSplitFile = 0 ; 
123 }
124
125 //____________________________________________________________________________ 
126 void AliPHOSSDigitizer::Init()
127 {
128   // Initialization: open root-file, allocate arrays for hits and sdigits,
129   // attach task SDigitizer to the list of PHOS tasks
130   // 
131   // Initialization can not be done in the default constructor
132   //============================================================= YS
133   //  The initialisation is now done by AliPHOSGetter
134   
135   if( strcmp(GetTitle(), "") == 0 )
136     SetTitle("galice.root") ;
137
138   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName()) ;     
139   if ( gime == 0 ) {
140     cerr << "ERROR: AliPHOSSDigitizer::Init -> Could not obtain the Getter object !" << endl ; 
141     return ;
142   } 
143   
144   gime->PostSDigits( GetName(), GetTitle() ) ; 
145
146   TString sdname(GetName() );
147   sdname.Append(":") ;
148   sdname.Append(GetTitle() ) ;
149   SetName(sdname) ;
150   gime->PostSDigitizer(this) ;
151 }
152
153 //____________________________________________________________________________ 
154 void AliPHOSSDigitizer::InitParameters()
155
156   fA             = 0;
157   fB             = 10000000.;
158   fPrimThreshold = 0.01 ;
159   fSDigitsInRun  = 0 ;
160   fSplitFile     = 0 ; 
161 }
162
163 //____________________________________________________________________________
164 void AliPHOSSDigitizer::Exec(Option_t *option) 
165
166   // Collects all hits in the same active volume into digit
167   if( strcmp(GetName(), "") == 0 )
168     Init() ;
169   
170   if (strstr(option, "print") ) {
171     Print("") ; 
172     return ; 
173   }
174
175   if(strstr(option,"tim"))
176     gBenchmark->Start("PHOSSDigitizer");
177
178   //Check, if this branch already exits
179   gAlice->GetEvent(0) ;
180
181   TString sdname(GetName()) ;
182   sdname.Remove(sdname.Index(GetTitle())-1) ;
183   
184   if(gAlice->TreeS() ) {
185     TObjArray * lob = static_cast<TObjArray*>(gAlice->TreeS()->GetListOfBranches()) ;
186     TIter next(lob) ; 
187     TBranch * branch = 0 ;  
188     Bool_t phosfound = kFALSE, sdigitizerfound = kFALSE ; 
189     
190     while ( (branch = (static_cast<TBranch*>(next()))) && (!phosfound || !sdigitizerfound) ) {
191       TString thisName( GetName() ) ; 
192       TString branchName( branch->GetTitle() ) ; 
193       branchName.Append(":") ; 
194       if ( (strcmp(branch->GetName(), "PHOS")==0) && thisName.BeginsWith(branchName) )  
195         phosfound = kTRUE ;
196       else if ( (strcmp(branch->GetName(), "AliPHOSSDigitizer")==0) &&  thisName.BeginsWith(branchName) )
197         sdigitizerfound = kTRUE ; 
198     }
199     
200     if ( phosfound || sdigitizerfound ) {
201       cerr << "WARNING: AliPHOSSDigitizer::Exec -> SDigits and/or SDigitizer branch with name " << GetName() 
202            << " already exits" << endl ;
203       return ; 
204     }   
205   }  
206
207   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
208
209   Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ; 
210   
211   Int_t ievent ;
212   for(ievent = 0; ievent < nevents; ievent++){
213     gime->Event(ievent,"H") ;
214     const TClonesArray * hits = gime->Hits() ;
215     TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
216     sdigits->Clear();
217     Int_t nSdigits = 0 ;
218     
219     
220     //Now make SDigits from hits, for PHOS it is the same, so just copy    
221
222     Int_t Nprim =  (Int_t)  (gAlice->TreeH())->GetEntries(); 
223     // Attention Nprim is the number of primaries tracked by Geant and this number could be different to the number of Primaries in TreeK;
224     Int_t iprim;
225     for (iprim=0; iprim<Nprim; iprim++) { 
226       //=========== Get the PHOS branch from Hits Tree for the Primary iprim
227       gime->Track(iprim) ;
228       Int_t i;
229       for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
230         AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
231         // Assign primary number only if contribution is significant
232         
233         if( hit->GetEnergy() > fPrimThreshold)
234           new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
235                                                   Digitize(hit->GetEnergy()), hit->GetTime()) ;
236         else
237           new((*sdigits)[nSdigits]) AliPHOSDigit( -1              , hit->GetId(), 
238                                                    Digitize(hit->GetEnergy()), hit->GetTime()) ;
239         nSdigits++ ;    
240         
241       }
242
243     } // loop over iprim
244     
245     sdigits->Sort() ;
246     
247     nSdigits = sdigits->GetEntriesFast() ;
248     fSDigitsInRun += nSdigits ;  
249     sdigits->Expand(nSdigits) ;
250
251     Int_t i ;
252     for (i = 0 ; i < nSdigits ; i++) { 
253       AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ; 
254       digit->SetIndexInList(i) ;     
255     }
256
257     if(gAlice->TreeS() == 0)
258       gAlice->MakeTree("S", fSplitFile);
259     
260     //First list of sdigits
261     Int_t bufferSize = 32000 ;    
262     TBranch * sdigitsBranch = gAlice->TreeS()->Branch("PHOS",&sdigits,bufferSize);
263     sdigitsBranch->SetTitle(sdname);
264       
265     //Next - SDigitizer
266     Int_t splitlevel = 0 ;
267     AliPHOSSDigitizer * sd = this ;
268     TBranch * sdigitizerBranch = gAlice->TreeS()->Branch("AliPHOSSDigitizer","AliPHOSSDigitizer",
269                                                &sd,bufferSize,splitlevel); 
270     sdigitizerBranch->SetTitle(sdname);
271   
272     sdigitsBranch->Fill() ; 
273     sdigitizerBranch->Fill() ;
274     gAlice->TreeS()->AutoSave() ;
275     
276     if(strstr(option,"deb"))
277       PrintSDigits(option) ;
278   
279   }
280
281   if(strstr(option,"tim")){
282     gBenchmark->Stop("PHOSSDigitizer");
283     cout << "AliPHOSSDigitizer:" << endl ;
284     cout << "   took " << gBenchmark->GetCpuTime("PHOSSDigitizer") << " seconds for SDigitizing " 
285          <<  gBenchmark->GetCpuTime("PHOSSDigitizer")/nevents << " seconds per event " << endl ;
286     cout << endl ;
287   }
288   
289   
290 }
291
292 //__________________________________________________________________
293 void AliPHOSSDigitizer::SetSDigitsBranch(const char * title )
294 {
295   // Setting title to branch SDigits 
296
297   TString stitle(title) ;
298
299   // check if branch with title already exists
300   TBranch * sdigitsBranch    = 
301     static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("PHOS")) ; 
302   TBranch * sdigitizerBranch =  
303     static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("AliPHOSSDigitizer")) ;
304   const char * sdigitsTitle    = sdigitsBranch ->GetTitle() ;  
305   const char * sdigitizerTitle = sdigitizerBranch ->GetTitle() ;
306   if ( stitle.CompareTo(sdigitsTitle)==0 || stitle.CompareTo(sdigitizerTitle)==0 ){
307     cerr << "ERROR: AliPHOSSdigitizer::SetSDigitsBranch -> Cannot overwrite existing branch with title " << title << endl ;
308     return ;
309   }
310   
311   cout << "AliPHOSSdigitizer::SetSDigitsBranch -> Changing SDigits file from " << GetName() << " to " << title << endl ;
312
313   SetName(title) ; 
314     
315   // Post to the WhiteBoard
316   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
317   gime->PostSDigits( title, GetTitle()) ; 
318 }
319
320 //__________________________________________________________________
321 void AliPHOSSDigitizer::SetSplitFile(const TString splitFileName) 
322 {
323   // Diverts the SDigits in a file separate from the hits file
324   
325   TDirectory * cwd = gDirectory ;
326   
327   if ( !(gAlice->GetTreeSFileName() == splitFileName) ) {
328     if (gAlice->GetTreeSFile() ) 
329       gAlice->GetTreeSFile()->Close() ; 
330   }
331   
332   fSplitFile = gAlice->InitTreeFile("S",splitFileName.Data());
333   fSplitFile->cd() ; 
334   gAlice->Write(0, TObject::kOverwrite);
335   
336   TTree *treeE  = gAlice->TreeE();
337   if (!treeE) {
338     cerr << "ERROR: AliPHOSSDigitizer::SetSPlitFile -> No TreeE found "<<endl;
339     abort() ;
340   }      
341   
342   // copy TreeE
343     AliHeader *header = new AliHeader();
344     treeE->SetBranchAddress("Header", &header);
345     treeE->SetBranchStatus("*",1);
346     TTree *treeENew =  treeE->CloneTree();
347     treeENew->Write(0, TObject::kOverwrite);
348   
349
350   // copy AliceGeom
351     TGeometry *AliceGeom = static_cast<TGeometry*>(cwd->Get("AliceGeom"));
352     if (!AliceGeom) {
353       cerr << "ERROR: AliPHOSSDigitizer::SetSPlitFile -> AliceGeom was not found in the input file "<<endl;
354       abort() ;
355     }
356     AliceGeom->Write(0, TObject::kOverwrite);
357
358   gAlice->MakeTree("S",fSplitFile);
359   cwd->cd() ; 
360   cout << "INFO: AliPHOSSDigitizer::SetSPlitMode -> SDigits will be stored in " << splitFileName.Data() << endl ; 
361
362 }
363
364 //__________________________________________________________________
365 void AliPHOSSDigitizer::Print(Option_t* option)const
366 {
367   // Prints parameters of SDigitizer
368   cout << "------------------- "<< GetName() << " -------------" << endl ;
369   cout << "   Writing SDigits to branch with title  " << GetName() << endl ;
370   cout << "   with digitization parameters  A = " << fA << endl ;
371   cout << "                                 B = " << fB << endl ;
372   cout << "   Threshold for Primary assignment= " << fPrimThreshold << endl ; 
373   cout << "---------------------------------------------------"<<endl ;
374   
375 }
376 //__________________________________________________________________
377 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
378 {
379   // Equal operator.
380   // SDititizers are equal if their pedestal, slope and threshold are equal
381
382   if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
383     return kTRUE ;
384   else
385     return kFALSE ;
386 }
387 //__________________________________________________________________
388 //__________________________________________________________________
389 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
390 {
391   // Prints list of digits produced in the current pass of AliPHOSDigitizer
392
393
394   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
395   TString sdname(GetName()) ;
396   sdname.Remove(sdname.Index(GetTitle())-1) ;
397   const TClonesArray * sdigits = gime->SDigits(sdname.Data()) ; 
398
399   cout << "AliPHOSSDigitiser: event " << gAlice->GetEvNumber() << endl ;
400   cout << "      Number of entries in SDigits list " << sdigits->GetEntriesFast() << endl ;
401   cout << endl ;
402   if(strstr(option,"all")||strstr(option,"EMC")){
403     
404     //loop over digits
405     AliPHOSDigit * digit;
406     cout << "EMC sdigits " << endl ;
407     cout << "Digit Id    Amplitude     Index     Nprim  Primaries list " <<  endl;      
408     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
409     Int_t index ;
410     for (index = 0 ; (index < sdigits->GetEntriesFast()) && 
411          (((AliPHOSDigit * )  sdigits->At(index))->GetId() <= maxEmc) ; index++) {
412       digit = (AliPHOSDigit * )  sdigits->At(index) ;
413       if(digit->GetNprimary() == 0) continue;
414       cout << setw(6)  <<  digit->GetId() << "   "  <<  setw(10)  <<  digit->GetAmp() <<   "    "  
415            << setw(6)  <<  digit->GetIndexInList() << "    "   
416            << setw(5)  <<  digit->GetNprimary() <<"    ";
417       
418       Int_t iprimary;
419       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
420         cout << setw(5)  <<  digit->GetPrimary(iprimary+1) << "    ";
421       cout << endl;      
422     }    
423     cout << endl;
424   }
425
426   if(strstr(option,"all")||strstr(option,"CPV")){
427     
428     //loop over CPV digits
429     AliPHOSDigit * digit;
430     cout << "CPV sdigits " << endl ;
431     cout << "Digit Id    Amplitude     Index     Nprim  Primaries list " <<  endl;      
432     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
433     Int_t index ;
434     for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
435       digit = (AliPHOSDigit * )  sdigits->At(index) ;
436       if(digit->GetId() > maxEmc){
437         cout << setw(6)  <<  digit->GetId() << "   "  <<        setw(10)  <<  digit->GetAmp() <<   "    "  
438              << setw(6)  <<  digit->GetIndexInList() << "    "   
439              << setw(5)  <<  digit->GetNprimary() <<"    ";
440         
441         Int_t iprimary;
442         for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
443           cout << setw(5)  <<  digit->GetPrimary(iprimary+1) << "    ";
444         cout << endl;    
445       }    
446     }
447   }
448
449 }
450
451 //____________________________________________________________________________ 
452 void AliPHOSSDigitizer::UseHitsFrom(const char * filename)
453 {
454   SetTitle(filename) ; 
455   Init() ; 
456 }