]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSSDigitizer.cxx
Transition to NewIO
[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 "TROOT.h"
56 #include "TBenchmark.h"
57
58 // --- Standard library ---
59
60 // --- AliRoot header files ---
61 #include "AliRun.h"
62 #include "AliPHOSDigit.h"
63 #include "AliPHOSGetter.h"
64 #include "AliPHOSHit.h"
65 #include "AliPHOSSDigitizer.h"
66
67 ClassImp(AliPHOSSDigitizer)
68
69            
70 //____________________________________________________________________________ 
71   AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("","")
72 {
73   // ctor
74   InitParameters() ;
75   fDefaultInit = kTRUE ; 
76 }
77
78 //____________________________________________________________________________ 
79 AliPHOSSDigitizer::AliPHOSSDigitizer(const char * alirunFileName, const char * eventFolderName):
80   TTask("PHOS"+AliConfig::fgkSDigitizerTaskName, alirunFileName),
81   fEventFolderName(eventFolderName)
82 {
83
84   // ctor
85   InitParameters() ; 
86   Init();
87   fDefaultInit = kFALSE ; 
88 }
89
90 //____________________________________________________________________________ 
91 AliPHOSSDigitizer::AliPHOSSDigitizer(const AliPHOSSDigitizer & sd) {
92   //cpy ctor 
93
94   fA             = sd.fA ;
95   fB             = sd.fB ;
96   fPrimThreshold = sd.fPrimThreshold ;
97   fSDigitsInRun  = sd.fSDigitsInRun ;
98   SetName(sd.GetName()) ; 
99   SetTitle(sd.GetTitle()) ; 
100   fEventFolderName = sd.fEventFolderName;
101 }
102
103
104 //____________________________________________________________________________ 
105 void AliPHOSSDigitizer::Init()
106 {
107   // Uses the getter to access the required files
108   
109   fInit = kTRUE ; 
110   
111   AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());  
112   if ( gime == 0 ) {
113     Fatal("Init" ,"Could not obtain the Getter object for file %s and event %s !", GetTitle(), fEventFolderName.Data()) ;  
114     return ;
115   } 
116   
117   TString opt("SDigits") ; 
118   if(gime->VersionExists(opt) ) { 
119     Error( "Init", "Give a version name different from %s", fEventFolderName.Data() ) ;
120     fInit = kFALSE ; 
121   }
122
123   gime->PostSDigitizer(this);
124   gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
125   
126 }
127
128 //____________________________________________________________________________ 
129 void AliPHOSSDigitizer::InitParameters()
130
131   // initializes the parameters for digitization
132   fA             = 0;
133   fB             = 10000000.;
134   fPrimThreshold = 0.01 ;
135   fSDigitsInRun  = 0 ;
136 }
137
138 //____________________________________________________________________________
139 void AliPHOSSDigitizer::Exec(Option_t *option) 
140
141   // Collects all hits in the same active volume into digit
142   
143   if (strstr(option, "print") ) {
144     Print() ; 
145     return ; 
146   }
147
148   if(strstr(option,"tim"))
149     gBenchmark->Start("PHOSSDigitizer");
150   
151   AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
152
153   //switch off reloading of this task while getting event
154   if (!fInit) { // to prevent overwrite existing file
155     Error( "Exec", "Give a version name different from %s", fEventFolderName.Data() ) ;
156     return ;
157   }
158
159
160   Int_t nevents = gime->MaxEvent() ; 
161   Int_t ievent ;
162   for(ievent = 0; ievent < nevents; ievent++){
163
164     gime->Event(ievent,"H") ;
165
166     TTree * treeS = gime->TreeS(); 
167     TClonesArray * hits = gime->Hits() ;
168     TClonesArray * sdigits = gime->SDigits() ;
169     sdigits->Clear();
170     Int_t nSdigits = 0 ;
171     //Now make SDigits from hits, for PHOS it is the same, so just copy    
172     Int_t nPrim =  static_cast<Int_t>((gime->TreeH())->GetEntries()) ; 
173     // Attention nPrim is the number of primaries tracked by Geant 
174     // and this number could be different to the number of Primaries in TreeK;
175     Int_t iprim ;
176
177     for (iprim = 0 ; iprim < nPrim ; iprim ++) { 
178       //=========== Get the PHOS branch from Hits Tree for the Primary iprim
179       gime->Track(iprim) ;
180       Int_t i;
181       for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
182         AliPHOSHit * hit = dynamic_cast<AliPHOSHit *>(hits->At(i)) ;
183         // Assign primary number only if contribution is significant
184         
185         if( hit->GetEnergy() > fPrimThreshold)
186           new((*sdigits)[nSdigits]) AliPHOSDigit(hit->GetPrimary(),hit->GetId(),
187                                                  Digitize(hit->GetEnergy()), hit->GetTime()) ;
188         else
189           new((*sdigits)[nSdigits]) AliPHOSDigit( -1              , hit->GetId(), 
190                                                   Digitize(hit->GetEnergy()), hit->GetTime()) ;
191         nSdigits++ ;    
192         
193       }
194  
195     } // loop over iprim
196
197     sdigits->Sort() ;
198
199     nSdigits = sdigits->GetEntriesFast() ;
200
201     fSDigitsInRun += nSdigits ;  
202     sdigits->Expand(nSdigits) ;
203
204     Int_t i ;
205     for (i = 0 ; i < nSdigits ; i++) { 
206       AliPHOSDigit * digit = dynamic_cast<AliPHOSDigit *>(sdigits->At(i)) ; 
207       digit->SetIndexInList(i) ;     
208     }
209
210     //Now write SDigits
211
212     
213     //First list of sdigits
214
215     Int_t bufferSize = 32000 ;
216     TBranch * sdigitsBranch = treeS->Branch("PHOS",&sdigits,bufferSize);
217
218     sdigitsBranch->Fill() ;
219
220     gime->WriteSDigits("OVERWRITE");
221
222     //Next - SDigitizer
223
224     gime->WriteSDigitizer("OVERWRITE");
225
226     if(strstr(option,"deb"))
227       PrintSDigits(option) ;
228   }
229   
230   Unload();
231
232   gime->PhosLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
233   
234   if(strstr(option,"tim")){
235     gBenchmark->Stop("PHOSSDigitizer");
236     Info("Exec","   took %f seconds for SDigitizing  %f seconds per event",
237          gBenchmark->GetCpuTime("PHOSSDigitizer"), gBenchmark->GetCpuTime("PHOSSDigitizer")/nevents) ;
238   }
239 }
240
241 //__________________________________________________________________
242 void AliPHOSSDigitizer::Print()const
243 {
244   // Prints parameters of SDigitizer
245   Info("Print", "\n------------------- %s -------------", GetName() ) ; 
246   printf("   Writing SDigits to branch with title  %s\n", fEventFolderName.Data()) ;
247   printf("   with digitization parameters  A = %f\n", fA) ; 
248   printf("                                 B = %f\n", fB) ;
249   printf("   Threshold for Primary assignment= %f\n", fPrimThreshold)  ; 
250   printf("---------------------------------------------------\n") ;
251   
252 }
253
254 //__________________________________________________________________
255 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const
256 {
257   // Equal operator.
258   // SDititizers are equal if their pedestal, slope and threshold are equal
259
260   if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
261     return kTRUE ;
262   else
263     return kFALSE ;
264 }
265
266 //__________________________________________________________________
267 void AliPHOSSDigitizer::PrintSDigits(Option_t * option)
268 {
269   // Prints list of digits produced in the current pass of AliPHOSDigitizer
270
271
272   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
273   const TClonesArray * sdigits = gime->SDigits() ;
274
275   TString message ; 
276   message  = "\nAliPHOSSDigitiser: event " ;
277   message += gAlice->GetEvNumber(); 
278   message += "\n      Number of entries in SDigits list " ;  
279   message += sdigits->GetEntriesFast() ; 
280   char * tempo = new char[8192]; 
281   
282   if(strstr(option,"all")||strstr(option,"EMC")){
283     
284     //loop over digits
285     AliPHOSDigit * digit;
286     message += "\nEMC sdigits\n" ;
287     message += "\n   Id  Amplitude    Time          Index Nprim: Primaries list \n" ;    
288     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
289     Int_t index ;
290     for (index = 0 ; (index < sdigits->GetEntriesFast()) && 
291          ((dynamic_cast<AliPHOSDigit *> (sdigits->At(index)))->GetId() <= maxEmc) ; index++) {
292       digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
293       if(digit->GetNprimary() == 0) 
294         continue;
295       sprintf(tempo, "\n%6d  %8d    %6.5e %4d      %2d :",
296               digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;  
297       message += tempo ; 
298       Int_t iprimary;
299       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
300         sprintf(tempo, "%d ",digit->GetPrimary(iprimary+1) ) ; 
301         message += tempo ; 
302       }  
303     }    
304   }
305
306   if(strstr(option,"all")||strstr(option,"CPV")){
307     
308     //loop over CPV digits
309     AliPHOSDigit * digit;
310     
311     message += "CPV sdigits\n" ;
312     message += "\n   Id  Amplitude    Index Nprim: Primaries list \n" ;    
313     Int_t maxEmc = gime->PHOSGeometry()->GetNModules()*gime->PHOSGeometry()->GetNCristalsInModule() ;
314     Int_t index ;
315     for (index = 0 ; index < sdigits->GetEntriesFast(); index++) {
316       digit = dynamic_cast<AliPHOSDigit *>( sdigits->At(index) ) ;
317       if(digit->GetId() > maxEmc){
318         sprintf(tempo, "\n%6d  %8d    %4d      %2d :",
319                 digit->GetId(), digit->GetAmp(), digit->GetIndexInList(), digit->GetNprimary()) ;  
320         message += tempo ; 
321         Int_t iprimary;
322         for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
323           sprintf(tempo, "%d ",digit->GetPrimary(iprimary+1) ) ; 
324           message += tempo ; 
325         }
326       }    
327     }
328   }
329   delete []tempo ; 
330   Info("PrintSDigits", message.Data() ) ;
331 }
332
333 //____________________________________________________________________________ 
334 void AliPHOSSDigitizer::Unload() const
335 {
336   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
337   AliPHOSLoader * loader = gime->PhosLoader() ; 
338   loader->UnloadHits() ; 
339   loader->UnloadSDigits() ; 
340 }