]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALSDigitizer.cxx
Macro to plot pathlengths of back-to-back jets. (A. Dainese)
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALSDigitizer.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 /* $Id$ */
17
18 //_________________________________________________________________________
19 // This is a TTask that makes SDigits out of Hits
20 // A Summable Digits is the sum of all hits originating 
21 // from one in one tower of the EMCAL 
22 // A threshold for assignment of the primary to SDigit is applied 
23 // SDigits are written to TreeS, branch "EMCAL"
24 // AliEMCALSDigitizer with all current parameters is written 
25 // to TreeS branch "AliEMCALSDigitizer".
26 // Both branches have the same title. If necessary one can produce 
27 // another set of SDigits with different parameters. Two versions
28 // can be distunguished using titles of the branches.
29 // User case:
30 // root [0] AliEMCALSDigitizer * s = new AliEMCALSDigitizer("galice.root")
31 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
32 // root [1] s->ExecuteTask()
33 //             // Makes SDigitis for all events stored in galice.root
34 // root [2] s->SetPedestalParameter(0.001)
35 //             // One can change parameters of digitization
36 // root [3] s->SetSDigitsBranch("Redestal 0.001")
37 //             // and write them into the new branch
38 // root [4] s->ExeciteTask("deb all tim")
39 //             // available parameters:
40 //             deb - print # of produced SDigitis
41 //             deb all  - print # and list of produced SDigits
42 //             tim - print benchmarking information
43 //
44 //*-- Author : Sahal Yacoob (LBL)
45 // based on  : AliPHOSSDigitzer 
46 // Modif: 
47 //  August 2002 Yves Schutz: clone PHOS as closely as possible and intoduction
48 //                           of new  IO (à la PHOS)
49 //////////////////////////////////////////////////////////////////////////////
50
51 // --- ROOT system ---
52 #include "TBenchmark.h"
53
54 // --- Standard library ---
55 #include "stdlib.h"
56
57 // --- AliRoot header files ---
58 #include "AliEMCALDigit.h"
59 #include "AliEMCALGetter.h"
60 #include "AliEMCALHit.h"
61 #include "AliEMCALSDigitizer.h"
62 #include "AliEMCALGeometry.h"
63
64 ClassImp(AliEMCALSDigitizer)
65            
66 //____________________________________________________________________________ 
67   AliEMCALSDigitizer::AliEMCALSDigitizer():TTask("","") 
68 {
69   // ctor
70   fFirstEvent = fLastEvent  = 0 ;  
71   fDefaultInit = kTRUE ; 
72 }
73
74 //____________________________________________________________________________ 
75 AliEMCALSDigitizer::AliEMCALSDigitizer(const char * alirunFileName, const char * eventFolderName):
76   TTask("EMCAL"+AliConfig::fgkSDigitizerTaskName, alirunFileName),
77   fEventFolderName(eventFolderName)
78 {
79   // ctor
80   fFirstEvent = fLastEvent  = 0 ; // runs one event by defaut  
81   Init();
82   InitParameters() ; 
83   fDefaultInit = kFALSE ; 
84 }
85
86
87 //____________________________________________________________________________ 
88 AliEMCALSDigitizer::AliEMCALSDigitizer(const AliEMCALSDigitizer & sd) : TTask(sd) {
89   //cpy ctor 
90
91   fFirstEvent    = sd.fFirstEvent ; 
92   fLastEvent     = sd.fLastEvent ;
93   fA             = sd.fA ;
94   fB             = sd.fB ;
95   fECPrimThreshold  = sd.fECPrimThreshold ;
96   fSDigitsInRun  = sd.fSDigitsInRun ;
97   SetName(sd.GetName()) ; 
98   SetTitle(sd.GetTitle()) ; 
99   fEventFolderName = sd.fEventFolderName;
100 }
101
102
103 //____________________________________________________________________________ 
104 AliEMCALSDigitizer::~AliEMCALSDigitizer() {
105   // dtor
106   AliEMCALGetter * gime = 
107     AliEMCALGetter::Instance(GetTitle(), fEventFolderName.Data());  
108   gime->EmcalLoader()->CleanSDigitizer();
109 }
110
111 //____________________________________________________________________________ 
112 void AliEMCALSDigitizer::Init(){
113   // Initialization: open root-file, allocate arrays for hits and sdigits,
114   // attach task SDigitizer to the list of EMCAL tasks
115   // 
116   // Initialization can not be done in the default constructor
117   //============================================================= YS
118   //  The initialisation is now done by the getter
119
120   fInit = kTRUE ; 
121    
122   AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName.Data());  
123   if ( gime == 0 ) {
124     Fatal("Init", "Could not obtain the Getter objectfor file %s and event %s !", GetTitle(), fEventFolderName.Data()) ;  
125     return ;
126   } 
127   
128   TString opt("SDigits") ; 
129   if(gime->VersionExists(opt) ) { 
130     Error( "Init", "Give a version name different from %s", fEventFolderName.Data() ) ;
131     fInit = kFALSE ; 
132   }
133   
134   gime->PostSDigitizer(this);
135   gime->EmcalLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
136   
137 }
138
139 //____________________________________________________________________________ 
140 void AliEMCALSDigitizer::InitParameters()
141 {
142   // Initializes parameters
143   fA                      = 0;
144   fB                      = 10000000.;
145
146   AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
147   const AliEMCALGeometry * geom = gime->EMCALGeometry() ; 
148   if (geom->GetSampling() == 0.) {
149     Fatal("InitParameters", "Sampling factor not set !") ; 
150   }
151 //   else
152 //     Info("InitParameters", "Sampling factor set to %f", geom->GetSampling()) ; 
153   
154   // this threshold corresponds approximately to 100 MeV
155   fECPrimThreshold     = 100E-3;
156 }
157
158 //____________________________________________________________________________
159 void AliEMCALSDigitizer::Exec(Option_t *option) 
160
161   // Collects all hit of the same tower into digits
162   if (strstr(option, "print") ) {
163     Print() ; 
164     return ; 
165   }
166   
167   if(strstr(option,"tim"))
168     gBenchmark->Start("EMCALSDigitizer");
169
170   AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle()) ;
171  
172   //switch off reloading of this task while getting event
173   if (!fInit) { // to prevent overwrite existing file
174     Error( "Exec", "Give a version name different from %s", fEventFolderName.Data() ) ;
175     return ;
176     }
177
178   if (fLastEvent == -1) 
179     fLastEvent = gime->MaxEvent() - 1 ;
180   else 
181     fLastEvent = TMath::Min(fFirstEvent, gime->MaxEvent()); // only ine event at the time
182   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
183
184   Int_t ievent ;   
185   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
186     gime->Event(ievent,"H") ;  
187     TTree * treeS = gime->TreeS(); 
188     TClonesArray * hits = gime->Hits() ; 
189     TClonesArray * sdigits = gime->SDigits() ;
190     sdigits->Clear();
191     Int_t nSdigits = 0 ;
192     //Now make SDigits from hits, for EMCAL it is the same, so just copy    
193     Int_t nPrim =  static_cast<Int_t>((gime->TreeH())->GetEntries()) ; 
194     // Attention nPrim is the number of primaries tracked by Geant 
195     // and this number could be different to the number of Primaries in TreeK;
196     Int_t iprim ;
197     for ( iprim = 0 ; iprim < nPrim ; iprim++ ) { 
198       //=========== Get the EMCAL branch from Hits Tree for the Primary iprim
199       gime->Track(iprim) ;
200       Int_t i;
201       for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
202         AliEMCALHit * hit = dynamic_cast<AliEMCALHit*>(hits->At(i)) ;
203         AliEMCALDigit * curSDigit = 0 ;
204         AliEMCALDigit * sdigit = 0 ;
205         Bool_t newsdigit = kTRUE; 
206
207         // Assign primary number only if deposited energy is significant
208         AliEMCALGeometry * geom = gime->EMCALGeometry() ; 
209         if( geom->IsInECA(hit->GetId()) ){
210           if( hit->GetEnergy() >  fECPrimThreshold )
211             curSDigit =  new AliEMCALDigit( hit->GetPrimary(),
212                                             hit->GetIparent(), hit->GetId(), 
213                                             Digitize(hit->GetEnergy()), hit->GetTime() ) ;
214           else
215             curSDigit =  new AliEMCALDigit( -1               , 
216                                             -1               ,
217                                             hit->GetId(), 
218                                             Digitize(hit->GetEnergy()), hit->GetTime() ) ;
219         }
220         else{
221           newsdigit = kFALSE;
222         }
223
224         Int_t check = 0 ;
225
226         for(check= 0; check < nSdigits ; check++) {
227           sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(check)) ;
228           if(curSDigit != 0){
229             if( sdigit->GetId() == curSDigit->GetId()) { // Are we in the same ECAL tower ?              
230             *sdigit = *sdigit + *curSDigit;
231             newsdigit = kFALSE;
232             }
233           }
234         }
235         if (newsdigit) {
236           new((*sdigits)[nSdigits])  AliEMCALDigit(*curSDigit);
237           nSdigits++ ;  
238         }
239         delete curSDigit ; 
240       }  // loop over all hits (hit = deposited energy/entering particle)
241     } // loop over iprim
242     sdigits->Sort() ;
243     
244     nSdigits = sdigits->GetEntriesFast() ;
245     fSDigitsInRun += nSdigits ;  
246     sdigits->Expand(nSdigits) ;
247
248     Int_t nPrimarymax = -1 ; 
249     Int_t i ;
250     for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) { 
251       AliEMCALDigit * sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(i)) ;
252       sdigit->SetIndexInList(i) ;
253     }
254     
255     for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) {   
256       if (((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) > nPrimarymax)
257         nPrimarymax = ((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) ;
258     }
259     
260     // Now write SDigits    
261     //First list of sdigits
262
263     Int_t bufferSize = 32000 ;    
264     TBranch * sdigitsBranch = treeS->Branch("EMCAL",&sdigits,bufferSize);
265  
266     sdigitsBranch->Fill() ;
267
268     gime->WriteSDigits("OVERWRITE");
269     
270     //NEXT - SDigitizer
271
272     gime->WriteSDigitizer("OVERWRITE");
273     
274     if(strstr(option,"deb"))
275       PrintSDigits(option) ;  
276   }
277
278   Unload();
279   
280   gime->EmcalLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
281   
282   if(strstr(option,"tim")){
283     gBenchmark->Stop("EMCALSDigitizer"); 
284     printf("Exec: took %f seconds for SDigitizing %f seconds per event", 
285          gBenchmark->GetCpuTime("EMCALSDigitizer"), gBenchmark->GetCpuTime("EMCALSDigitizer")/nEvents ) ; 
286   }
287 }
288
289
290 //__________________________________________________________________
291 void AliEMCALSDigitizer::Print()const
292
293   // Prints parameters of SDigitizer
294   printf("Print: \n------------------- %s -------------", GetName() ) ; 
295   printf("   Writing SDigits to branch with title  %s\n", fEventFolderName.Data()) ;
296   printf("   with digitization parameters  A = %f\n", fA) ; 
297   printf("                                 B = %f\n", fB) ;
298   printf("   Threshold for EC Primary assignment= %f\n", fECPrimThreshold)  ; 
299   printf("---------------------------------------------------\n") ;
300
301 }
302
303 //__________________________________________________________________
304 Bool_t AliEMCALSDigitizer::operator==( AliEMCALSDigitizer const &sd )const
305 {
306   // Equal operator.
307   // SDititizers are equal if their pedestal, slope and threshold are equal
308   if( (fA==sd.fA)&&(fB==sd.fB)&&
309       (fECPrimThreshold==sd.fECPrimThreshold))
310     return kTRUE ;
311   else
312     return kFALSE ;
313 }
314
315 //__________________________________________________________________
316 void AliEMCALSDigitizer::PrintSDigits(Option_t * option)
317 {
318   //Prints list of digits produced at the current pass of AliEMCALDigitizer
319   
320   
321   AliEMCALGetter * gime = AliEMCALGetter::Instance() ; 
322   const TClonesArray * sdigits = gime->SDigits() ; 
323   
324   printf("\n") ;  
325   printf("event %i", gAlice->GetEvNumber()) ;
326   printf("\n      Number of entries in SDigits list %i", sdigits->GetEntriesFast()); 
327   if(strstr(option,"all")||strstr(option,"EMC")){
328     
329     //loop over digits
330     AliEMCALDigit * digit;
331     printf("\n   Id  Amplitude    Time          Index Nprim: Primaries list \n") ;    
332     Int_t index ;
333     char * tempo = new char[8192]; 
334     for (index = 0 ; index < sdigits->GetEntries() ; index++) {
335       digit = dynamic_cast<AliEMCALDigit *>( sdigits->At(index) ) ;
336       sprintf(tempo, "\n%6d  %8d    %6.5e %4d      %2d :",
337               digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;  
338       printf(tempo); 
339       
340       Int_t iprimary;
341       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
342         sprintf(tempo, "%d ",digit->GetPrimary(iprimary+1) ) ; 
343         printf(tempo); 
344       }          
345     }
346     delete tempo ;
347   }
348 }
349
350 //____________________________________________________________________________ 
351 void AliEMCALSDigitizer::Unload() const
352 {
353   // Unload Hits and SDigits from the folder
354   AliEMCALGetter * gime = AliEMCALGetter::Instance() ; 
355   AliEMCALLoader * loader = gime->EmcalLoader() ; 
356   loader->UnloadHits() ; 
357   loader->UnloadSDigits() ; 
358 }