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