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