5cd484e927f58955cd3600625acdde887de92706
[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 #include "TGeometry.h"
61
62 // --- Standard library ---
63 #include <iomanip.h>
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
75 ClassImp(AliPHOSSDigitizer)
76
77            
78 //____________________________________________________________________________ 
79   AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("","") 
80 {
81   // ctor
82   fA             = 0;
83   fB             = 10000000.;
84   fPrimThreshold = 0.01 ;
85   fSDigitsInRun  = 0 ; 
86 }
87
88 //____________________________________________________________________________ 
89 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * headerFile, const char * sDigitsTitle):TTask(sDigitsTitle, headerFile)
90 {
91   // ctor
92   fA             = 0;
93   fB             = 10000000.;
94   fPrimThreshold = 0.01 ;
95   fSDigitsInRun  = 0 ;
96   Init();
97 }
98
99
100 //____________________________________________________________________________ 
101 void AliPHOSSDigitizer::Init()
102 {
103   // Initialization: open root-file, allocate arrays for hits and sdigits,
104   // attach task SDigitizer to the list of PHOS tasks
105   // 
106   // Initialization can not be done in the default constructor
107   //============================================================= YS
108   //  The initialisation is now done by AliPHOSGetter
109   
110   if( strcmp(GetTitle(), "") == 0 )
111     SetTitle("galice.root") ;
112   
113   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), GetName()) ;     
114   if ( gime == 0 ) {
115     cerr << "ERROR: AliPHOSSDigitizer::Init -> Could not obtain the Getter object !" << endl ; 
116     return ;
117   } 
118   
119   gime->PostSDigits( GetName(), GetTitle() ) ; 
120
121   TString sdname(GetName() );
122   sdname.Append(":") ;
123   sdname.Append(GetTitle() ) ;
124   SetName(sdname) ;
125   gime->PostSDigitizer(this) ;
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
146   TString sdname(GetName()) ;
147   sdname.Remove(sdname.Index(GetTitle())-1) ;
148   
149   if(gAlice->TreeS() ) {
150     TObjArray * lob = static_cast<TObjArray*>(gAlice->TreeS()->GetListOfBranches()) ;
151     TIter next(lob) ; 
152     TBranch * branch = 0 ;  
153     Bool_t phosfound = kFALSE, sdigitizerfound = kFALSE ; 
154     
155     while ( (branch = (static_cast<TBranch*>(next()))) && (!phosfound || !sdigitizerfound) ) {
156       if ( (strcmp(branch->GetName(), "PHOS")==0) && (strcmp(branch->GetTitle(), GetName())==0) ) 
157         phosfound = kTRUE ;
158       
159       else if ( (strcmp(branch->GetName(), "AliPHOSSDigitizer")==0) && (strcmp(branch->GetTitle(), GetName())==0) ) 
160         sdigitizerfound = kTRUE ; 
161     }
162     
163     if ( phosfound || sdigitizerfound ) {
164       cerr << "WARNING: AliPHOSSDigitizer::Exec -> SDigits and/or SDigitizer branch with name " << GetName() 
165            << " already exits" << endl ;
166       return ; 
167     }   
168   }  
169
170   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
171
172   Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ; 
173   
174   Int_t ievent ;
175   for(ievent = 0; ievent < nevents; ievent++){
176     gime->Event(ievent,"H") ;
177     const TClonesArray * hits = gime->Hits() ;
178     TClonesArray * sdigits = gime->SDigits(sdname.Data()) ;
179     sdigits->Clear();
180     Int_t nSdigits = 0 ;
181     
182     
183     //Now make SDigits from hits, for PHOS it is the same, so just copy    
184     Int_t Nprim =  (Int_t)  (gAlice->TreeH())->GetEntries(); 
185     // Attention Nprim is the number of primaries tracked by Geant and this number could be different to the number of Primaries in TreeK;
186     Int_t iprim;
187     for (iprim=0; iprim<Nprim; iprim++) { 
188       //=========== Get the PHOS branch from Hits Tree for the Primary iprim
189       gime->Track(iprim) ;
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 iprim
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     //First list of sdigits
221     Int_t bufferSize = 32000 ;    
222     TBranch * sdigitsBranch = gAlice->TreeS()->Branch("PHOS",&sdigits,bufferSize);
223     sdigitsBranch->SetTitle(sdname);
224       
225     //Next - SDigitizer
226     Int_t splitlevel = 0 ;
227     AliPHOSSDigitizer * sd = this ;
228     TBranch * sdigitizerBranch = gAlice->TreeS()->Branch("AliPHOSSDigitizer","AliPHOSSDigitizer",
229                                                &sd,bufferSize,splitlevel); 
230     sdigitizerBranch->SetTitle(sdname);
231     sdigitsBranch->Fill() ; 
232
233     gAlice->TreeS()->AutoSave() ;
234     
235     if(strstr(option,"deb"))
236       PrintSDigits(option) ;
237     
238   }
239   
240   if(strstr(option,"tim")){
241     gBenchmark->Stop("PHOSSDigitizer");
242     cout << "AliPHOSSDigitizer:" << endl ;
243     cout << "   took " << gBenchmark->GetCpuTime("PHOSSDigitizer") << " seconds for SDigitizing " 
244          <<  gBenchmark->GetCpuTime("PHOSSDigitizer")/nevents << " seconds per event " << endl ;
245     cout << endl ;
246   }
247   
248   
249 }
250
251 //__________________________________________________________________
252 void AliPHOSSDigitizer::SetSDigitsBranch(const char * title )
253 {
254   // Setting title to branch SDigits 
255
256   TString stitle(title) ;
257
258   // check if branch with title already exists
259   TBranch * sdigitsBranch    = 
260     static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("PHOS")) ; 
261   TBranch * sdigitizerBranch =  
262     static_cast<TBranch*>(gAlice->TreeS()->GetListOfBranches()->FindObject("AliPHOSSDigitizer")) ;
263   const char * sdigitsTitle    = sdigitsBranch ->GetTitle() ;  
264   const char * sdigitizerTitle = sdigitizerBranch ->GetTitle() ;
265   if ( stitle.CompareTo(sdigitsTitle)==0 || stitle.CompareTo(sdigitizerTitle)==0 ){
266     cerr << "ERROR: AliPHOSSdigitizer::SetSDigitsBranch -> Cannot overwrite existing branch with title " << title << endl ;
267     return ;
268   }
269   
270   cout << "AliPHOSSdigitizer::SetSDigitsBranch -> Changing SDigits file from " << GetName() << " to " << title << endl ;
271
272   SetName(title) ; 
273     
274   // Post to the WhiteBoard
275   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
276   gime->PostSDigits( title, GetTitle()) ; 
277 }
278
279 //__________________________________________________________________
280 void AliPHOSSDigitizer::SetSplitFile(const TString splitFileName) const
281 {
282   // Diverts the SDigits in a file separate from the hits file
283   
284   TDirectory * cwd = gDirectory ;
285   TFile * splitFile = gAlice->InitTreeFile("S",splitFileName.Data());
286   splitFile->cd() ; 
287   gAlice->Write();
288   
289   TTree *treeE  = gAlice->TreeE();
290   if (!treeE) {
291     cerr<<"No TreeE found "<<endl;
292     abort() ;
293   }      
294   
295   // copy TreeE
296   AliHeader *header = new AliHeader();
297   treeE->SetBranchAddress("Header", &header);
298   treeE->SetBranchStatus("*",1);
299   TTree *treeENew =  treeE->CloneTree();
300   treeENew->Write();
301   
302   // copy AliceGeom
303   TGeometry *AliceGeom = static_cast<TGeometry*>(cwd->Get("AliceGeom"));
304   if (!AliceGeom) {
305     cerr<<"AliceGeom was not found in the input file "<<endl;
306     abort() ;
307   }
308   AliceGeom->Write();
309   cwd->cd() ; 
310   gAlice->MakeTree("S",splitFile);
311   cout << "INFO: AliPHOSSDigitizer::SetSPlitMode -> SDigits will be stored in " << splitFileName.Data() << endl ; 
312 }
313
314 //__________________________________________________________________
315 void AliPHOSSDigitizer::Print(Option_t* option)const
316 {
317   // Prints parameters of SDigitizer
318   cout << "------------------- "<< GetName() << " -------------" << endl ;
319   cout << "   Writing SDigits to branch with title  " << GetName() << endl ;
320   cout << "   with digitization parameters  A = " << fA << endl ;
321   cout << "                                 B = " << fB << endl ;
322   cout << "   Threshold for Primary assignment= " << fPrimThreshold << endl ; 
323   cout << "---------------------------------------------------"<<endl ;
324   
325 }
326 //__________________________________________________________________
327 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
328 {
329   // Equal operator.
330   // SDititizers are equal if their pedestal, slope and threshold are equal
331
332   if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
333     return kTRUE ;
334   else
335     return kFALSE ;
336 }
337 //__________________________________________________________________
338 //__________________________________________________________________
339 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
340 {
341   // Prints list of digits produced in the current pass of AliPHOSDigitizer
342
343
344   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
345   TString sdname(GetName()) ;
346   sdname.Remove(sdname.Index(GetTitle())-1) ;
347   const TClonesArray * sdigits = gime->SDigits(sdname.Data()) ; 
348
349   cout << "AliPHOSSDigitiser: event " << gAlice->GetEvNumber() << endl ;
350   cout << "      Number of entries in SDigits list " << sdigits->GetEntriesFast() << endl ;
351   cout << endl ;
352   if(strstr(option,"all")||strstr(option,"EMC")){
353     
354     //loop over digits
355     AliPHOSDigit * digit;
356     cout << "EMC sdigits " << endl ;
357     cout << "Digit Id    Amplitude     Index     Nprim  Primaries list " <<  endl;      
358     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
359     Int_t index ;
360     for (index = 0 ; (index < sdigits->GetEntriesFast()) && 
361          (((AliPHOSDigit * )  sdigits->At(index))->GetId() <= maxEmc) ; index++) {
362       digit = (AliPHOSDigit * )  sdigits->At(index) ;
363       if(digit->GetNprimary() == 0) continue;
364       cout << setw(6)  <<  digit->GetId() << "   "  <<  setw(10)  <<  digit->GetAmp() <<   "    "  
365            << setw(6)  <<  digit->GetIndexInList() << "    "   
366            << setw(5)  <<  digit->GetNprimary() <<"    ";
367       
368       Int_t iprimary;
369       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
370         cout << setw(5)  <<  digit->GetPrimary(iprimary+1) << "    ";
371       cout << endl;      
372     }    
373     cout << endl;
374   }
375
376   if(strstr(option,"all")||strstr(option,"CPV")){
377     
378     //loop over CPV digits
379     AliPHOSDigit * digit;
380     cout << "CPV sdigits " << endl ;
381     cout << "Digit Id    Amplitude     Index     Nprim  Primaries list " <<  endl;      
382     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
383     Int_t index ;
384     for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
385       digit = (AliPHOSDigit * )  sdigits->At(index) ;
386       if(digit->GetId() > maxEmc){
387         cout << setw(6)  <<  digit->GetId() << "   "  <<        setw(10)  <<  digit->GetAmp() <<   "    "  
388              << setw(6)  <<  digit->GetIndexInList() << "    "   
389              << setw(5)  <<  digit->GetNprimary() <<"    ";
390         
391         Int_t iprimary;
392         for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++)
393           cout << setw(5)  <<  digit->GetPrimary(iprimary+1) << "    ";
394         cout << endl;    
395       }    
396     }
397   }
398
399 }
400
401 //____________________________________________________________________________ 
402 void AliPHOSSDigitizer::UseHitsFrom(const char * filename)
403 {
404   SetTitle(filename) ; 
405   Init() ; 
406 }