]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALSDigitizer.cxx
final geometry and cleanup
[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   InitParameters() ; 
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   Init();
81   InitParameters() ; 
82   fDefaultInit = kFALSE ; 
83 }
84
85
86 //____________________________________________________________________________ 
87 AliEMCALSDigitizer::AliEMCALSDigitizer(const AliEMCALSDigitizer & sd) : TTask(sd) {
88   //cpy ctor 
89
90   fA             = sd.fA ;
91   fB             = sd.fB ;
92   fECPrimThreshold  = sd.fECPrimThreshold ;
93   fSDigitsInRun  = sd.fSDigitsInRun ;
94   SetName(sd.GetName()) ; 
95   SetTitle(sd.GetTitle()) ; 
96   fEventFolderName = sd.fEventFolderName;
97 }
98
99
100 //____________________________________________________________________________ 
101 void AliEMCALSDigitizer::Init(){
102   // Initialization: open root-file, allocate arrays for hits and sdigits,
103   // attach task SDigitizer to the list of EMCAL tasks
104   // 
105   // Initialization can not be done in the default constructor
106   //============================================================= YS
107   //  The initialisation is now done by the getter
108
109   fInit = kTRUE ; 
110    
111   AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName.Data());  
112   if ( gime == 0 ) {
113     Fatal("Init", "Could not obtain the Getter objectfor 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->EmcalLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
125   
126 }
127
128 //____________________________________________________________________________ 
129 void AliEMCALSDigitizer::InitParameters()
130 {
131   // Initializes parameters
132   fA                      = 0;
133   fB                      = 10000000.;
134
135   AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
136   const AliEMCALGeometry * geom = gime->EMCALGeometry() ; 
137   if (geom->GetSampling() == 0.) {
138     printf("InitParameters: Sampling factor not set !") ; 
139     abort() ;
140   }
141   else
142     printf("InitParameters: Sampling factor set to %f\n", geom->GetSampling()) ; 
143   
144   // this threshold corresponds approximately to 100 MeV
145   fECPrimThreshold     = 100E-3 / ( geom->GetSampling() * geom->GetNECLayers()) * geom->GetNECLayers() ;
146 }
147
148 //____________________________________________________________________________
149 void AliEMCALSDigitizer::Exec(Option_t *option) 
150
151   // Collects all hit of the same tower into digits
152   if (strstr(option, "print") ) {
153     Print() ; 
154     return ; 
155   }
156   
157   if(strstr(option,"tim"))
158     gBenchmark->Start("EMCALSDigitizer");
159
160   AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
161  
162   //switch off reloading of this task while getting event
163   if (!fInit) { // to prevent overwrite existing file
164     Error( "Exec", "Give a version name different from %s", fEventFolderName.Data() ) ;
165     return ;
166     }
167
168   Int_t nevents = gime->MaxEvent() ; 
169   Int_t ievent ;   
170   for(ievent = 0; ievent < nevents; ievent++){     
171     gime->Event(ievent,"H") ;  
172
173     TTree * treeS = gime->TreeS(); 
174     TClonesArray * hits = gime->Hits() ; 
175     TClonesArray * sdigits = gime->SDigits() ;
176     sdigits->Clear();
177     Int_t nSdigits = 0 ;
178     //Now make SDigits from hits, for EMCAL it is the same, so just copy    
179     Int_t nPrim =  static_cast<Int_t>((gime->TreeH())->GetEntries()) ; 
180     // Attention nPrim is the number of primaries tracked by Geant 
181     // and this number could be different to the number of Primaries in TreeK;
182     Int_t iprim ;
183     for ( iprim = 0 ; iprim < nPrim ; iprim++ ) { 
184       //=========== Get the EMCAL branch from Hits Tree for the Primary iprim
185       gime->Track(iprim) ;
186       Int_t i;
187       for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
188         AliEMCALHit * hit = dynamic_cast<AliEMCALHit*>(hits->At(i)) ;
189         AliEMCALDigit * curSDigit = 0 ;
190         AliEMCALDigit * sdigit = 0 ;
191         Bool_t newsdigit = kTRUE; 
192
193         // Assign primary number only if deposited energy is significant
194         AliEMCALGeometry * geom = gime->EMCALGeometry() ; 
195         if( geom->IsInECA(hit->GetId()) ){
196           if( hit->GetEnergy() >  fECPrimThreshold )
197             curSDigit =  new AliEMCALDigit( hit->GetPrimary(),
198                                             hit->GetIparent(), hit->GetId(), 
199                                             Digitize(hit->GetEnergy()), hit->GetTime() ) ;
200           else
201             curSDigit =  new AliEMCALDigit( -1               , 
202                                             -1               ,
203                                             hit->GetId(), 
204                                             Digitize(hit->GetEnergy()), hit->GetTime() ) ;
205         }
206         else{
207           newsdigit = kFALSE;
208         }
209
210         Int_t check = 0 ;
211
212         for(check= 0; check < nSdigits ; check++) {
213           sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(check)) ;
214           if(curSDigit != 0){
215             if( sdigit->GetId() == curSDigit->GetId()) { // Are we in the same ECAL tower ?              
216             *sdigit = *sdigit + *curSDigit;
217             newsdigit = kFALSE;
218             }
219           }
220         }
221         if (newsdigit) {
222           new((*sdigits)[nSdigits])  AliEMCALDigit(*curSDigit);
223           nSdigits++ ;  
224         }
225         delete curSDigit ; 
226       }  // loop over all hits (hit = deposited energy/entering particle)
227     } // loop over iprim
228     sdigits->Sort() ;
229     
230     nSdigits = sdigits->GetEntriesFast() ;
231     fSDigitsInRun += nSdigits ;  
232     sdigits->Expand(nSdigits) ;
233
234     Int_t nPrimarymax = -1 ; 
235     Int_t i ;
236     for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) { 
237       AliEMCALDigit * sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(i)) ;
238       sdigit->SetIndexInList(i) ;
239     }
240     
241     for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) {   
242       if (((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) > nPrimarymax)
243         nPrimarymax = ((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) ;
244     }
245     
246     // Now write SDigits    
247     //First list of sdigits
248
249     Int_t bufferSize = 32000 ;    
250     TBranch * sdigitsBranch = treeS->Branch("EMCAL",&sdigits,bufferSize);
251  
252     sdigitsBranch->Fill() ;
253
254     gime->WriteSDigits("OVERWRITE");
255     
256     //NEXT - SDigitizer
257
258     gime->WriteSDigitizer("OVERWRITE");
259     
260     if(strstr(option,"deb"))
261       PrintSDigits(option) ;  
262   }
263
264   Unload();
265   
266   gime->EmcalLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
267   
268   if(strstr(option,"tim")){
269     gBenchmark->Stop("EMCALSDigitizer"); 
270     printf("Exec: took %f seconds for SDigitizing %f seconds per event", 
271          gBenchmark->GetCpuTime("EMCALSDigitizer"), gBenchmark->GetCpuTime("EMCALSDigitizer") ) ; 
272   }
273 }
274
275
276 //__________________________________________________________________
277 void AliEMCALSDigitizer::Print()const
278
279   // Prints parameters of SDigitizer
280   printf("Print: \n------------------- %s -------------", GetName() ) ; 
281   printf("   Writing SDigits to branch with title  %s\n", fEventFolderName.Data()) ;
282   printf("   with digitization parameters  A = %f\n", fA) ; 
283   printf("                                 B = %f\n", fB) ;
284   printf("   Threshold for EC Primary assignment= %f\n", fECPrimThreshold)  ; 
285   printf("---------------------------------------------------\n") ;
286
287 }
288
289 //__________________________________________________________________
290 Bool_t AliEMCALSDigitizer::operator==( AliEMCALSDigitizer const &sd )const
291 {
292   // Equal operator.
293   // SDititizers are equal if their pedestal, slope and threshold are equal
294   if( (fA==sd.fA)&&(fB==sd.fB)&&
295       (fECPrimThreshold==sd.fECPrimThreshold))
296     return kTRUE ;
297   else
298     return kFALSE ;
299 }
300
301 //__________________________________________________________________
302 void AliEMCALSDigitizer::PrintSDigits(Option_t * option)
303 {
304   //Prints list of digits produced at the current pass of AliEMCALDigitizer
305   
306   
307   AliEMCALGetter * gime = AliEMCALGetter::Instance() ; 
308   const TClonesArray * sdigits = gime->SDigits() ; 
309   
310   printf("\n") ;  
311   printf("event %i", gAlice->GetEvNumber()) ;
312   printf("\n      Number of entries in SDigits list %i", sdigits->GetEntriesFast()); 
313   if(strstr(option,"all")||strstr(option,"EMC")){
314     
315     //loop over digits
316     AliEMCALDigit * digit;
317     printf("\n   Id  Amplitude    Time          Index Nprim: Primaries list \n") ;    
318     Int_t index ;
319     char * tempo = new char[8192]; 
320     for (index = 0 ; index < sdigits->GetEntries() ; index++) {
321       digit = dynamic_cast<AliEMCALDigit *>( sdigits->At(index) ) ;
322       sprintf(tempo, "\n%6d  %8d    %6.5e %4d      %2d :",
323               digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;  
324       printf(tempo); 
325       
326       Int_t iprimary;
327       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
328         sprintf(tempo, "%d ",digit->GetPrimary(iprimary+1) ) ; 
329         printf(tempo); 
330       }          
331     }
332     delete tempo ;
333   }
334 }
335
336 //____________________________________________________________________________ 
337 void AliEMCALSDigitizer::Unload() const
338 {
339   // Unload Hits and SDigits from the folder
340   AliEMCALGetter * gime = AliEMCALGetter::Instance() ; 
341   AliEMCALLoader * loader = gime->EmcalLoader() ; 
342   loader->UnloadHits() ; 
343   loader->UnloadSDigits() ; 
344 }