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