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