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