]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALSDigitizer.cxx
Small changes
[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 "TFile.h"
53 #include "TROOT.h"
54 #include "TBenchmark.h"
55
56 // --- Standard library ---
57
58 // --- AliRoot header files ---
59 #include "AliRun.h"
60 #include "AliEMCALDigit.h"
61 #include "AliEMCALGetter.h"
62 #include "AliEMCALHit.h"
63 #include "AliEMCALSDigitizer.h"
64 #include "AliEMCALGeometry.h"
65
66 ClassImp(AliEMCALSDigitizer)
67            
68 //____________________________________________________________________________ 
69   AliEMCALSDigitizer::AliEMCALSDigitizer():TTask("","") 
70 {
71   // ctor
72   InitParameters() ; 
73   fDefaultInit = kTRUE ; 
74 }
75
76 //____________________________________________________________________________ 
77 AliEMCALSDigitizer::AliEMCALSDigitizer(const char * alirunFileName, const char * eventFolderName):
78   TTask("EMCAL"+AliConfig::fgkSDigitizerTaskName, alirunFileName),
79   fEventFolderName(eventFolderName)
80 {
81   // ctor
82   Init();
83   InitParameters() ; 
84   fDefaultInit = kFALSE ; 
85 }
86
87
88 //____________________________________________________________________________ 
89 AliEMCALSDigitizer::AliEMCALSDigitizer(const AliEMCALSDigitizer & sd) : TTask(sd) {
90   //cpy ctor 
91
92   fA             = sd.fA ;
93   fB             = sd.fB ;
94   fPREPrimThreshold = sd.fPREPrimThreshold ;
95   fECPrimThreshold  = sd.fECPrimThreshold ;
96   fHCPrimThreshold  = sd.fHCPrimThreshold ;
97   fSDigitsInRun  = sd.fSDigitsInRun ;
98   SetName(sd.GetName()) ; 
99   SetTitle(sd.GetTitle()) ; 
100   fEventFolderName = sd.fEventFolderName;
101 }
102
103
104 //____________________________________________________________________________ 
105 void AliEMCALSDigitizer::Init(){
106   // Initialization: open root-file, allocate arrays for hits and sdigits,
107   // attach task SDigitizer to the list of EMCAL tasks
108   // 
109   // Initialization can not be done in the default constructor
110   //============================================================= YS
111   //  The initialisation is now done by the getter
112
113   fInit = kTRUE ; 
114    
115   AliEMCALGetter * gime = AliEMCALGetter::Instance(GetTitle(), fEventFolderName.Data());  
116   if ( gime == 0 ) {
117     Fatal("Init", "Could not obtain the Getter objectfor file %s and event %s !", GetTitle(), fEventFolderName.Data()) ;  
118     return ;
119   } 
120   
121   TString opt("SDigits") ; 
122   if(gime->VersionExists(opt) ) { 
123     Error( "Init", "Give a version name different from %s", fEventFolderName.Data() ) ;
124     fInit = kFALSE ; 
125   }
126   
127   gime->PostSDigitizer(this);
128   gime->EmcalLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
129   
130 }
131
132 //____________________________________________________________________________ 
133 void AliEMCALSDigitizer::InitParameters()
134 {
135   fA                      = 0;
136   fB                      = 10000000.;
137
138   AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
139   const AliEMCALGeometry * geom = gime->EMCALGeometry() ; 
140   if (geom->GetSampling() == 0.) {
141     Error("InitParameters", "Sampling factor not set !") ; 
142     abort() ;
143   }
144   else
145     Info("InitParameters", "Sampling factor set to %f\n", geom->GetSampling()) ; 
146   
147   // this threshold corresponds approximately to 100 MeV
148   fECPrimThreshold     = 100E-3 / ( geom->GetSampling() * ( geom->GetNPRLayers() + geom->GetNECLayers()) ) * geom->GetNECLayers() ;
149   fPREPrimThreshold    = 100E-3 / ( geom->GetSampling() * ( geom->GetNPRLayers() + geom->GetNECLayers()) ) * geom->GetNPRLayers() ; 
150   fHCPrimThreshold     = fECPrimThreshold/5. ; // 5 is totally arbitrary
151
152 }
153
154 //____________________________________________________________________________
155 void AliEMCALSDigitizer::Exec(Option_t *option) 
156
157   // Collects all hits in the section (PRE/ECAL/HCAL) of the same tower into digit
158   
159   if (strstr(option, "print") ) {
160     Print("") ; 
161     return ; 
162   }
163   
164   if(strstr(option,"tim"))
165     gBenchmark->Start("EMCALSDigitizer");
166
167   AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
168  
169   //switch off reloading of this task while getting event
170   if (!fInit) { // to prevent overwrite existing file
171     Error( "Exec", "Give a version name different from %s", fEventFolderName.Data() ) ;
172     return ;
173   }
174   
175   Int_t nevents = gime->MaxEvent() ; 
176   Int_t ievent ;   
177   for(ievent = 0; ievent < nevents; ievent++){     
178   
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
192     for ( iprim = 0 ; iprim < nPrim ; iprim++ ) { 
193       //=========== Get the EMCAL branch from Hits Tree for the Primary iprim
194       gime->Track(iprim) ;
195       Int_t i;
196       for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
197         AliEMCALHit * hit = dynamic_cast<AliEMCALHit*>(hits->At(i)) ;
198         AliEMCALDigit * curSDigit = 0 ;
199         AliEMCALDigit * sdigit = 0 ;
200         Bool_t newsdigit = kTRUE; 
201
202         // Assign primary number only if deposited energy is significant
203
204         AliEMCALGeometry * geom = gime->EMCALGeometry() ; 
205
206         if( geom->IsInPRE(hit->GetId()) )  
207           if( hit->GetEnergy() > fPREPrimThreshold )
208             curSDigit =  new AliEMCALDigit( hit->GetPrimary(),
209                                             hit->GetIparent(), hit->GetId(), 
210                                             Digitize(hit->GetEnergy()), hit->GetTime() ) ;
211           else
212             curSDigit =  new AliEMCALDigit( -1               , 
213                                             -1               ,
214                                             hit->GetId(), 
215                                             Digitize(hit->GetEnergy()), hit->GetTime() ) ;
216         else if( geom->IsInECA(hit->GetId()) )
217           if( hit->GetEnergy() >  fECPrimThreshold )
218             curSDigit =  new AliEMCALDigit( hit->GetPrimary(),
219                                             hit->GetIparent(), hit->GetId(), 
220                                             Digitize(hit->GetEnergy()), hit->GetTime() ) ;
221           else
222             curSDigit =  new AliEMCALDigit( -1               , 
223                                             -1               ,
224                                             hit->GetId(), 
225                                             Digitize(hit->GetEnergy()), hit->GetTime() ) ;
226         else if( geom->IsInHCA(hit->GetId()) )
227           if( hit->GetEnergy() >  fHCPrimThreshold )
228             
229             curSDigit =  new AliEMCALDigit( hit->GetPrimary(),
230                                             hit->GetIparent(), hit->GetId(), 
231                                             Digitize(hit->GetEnergy()), hit->GetTime() ) ;
232           else
233             curSDigit =  new AliEMCALDigit( -1               , 
234                                             -1               ,
235                                             hit->GetId(), 
236                                             Digitize(hit->GetEnergy()), hit->GetTime() ) ;
237         
238         Int_t check = 0 ;
239         for(check= 0; check < nSdigits ; check++) {
240           sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(check)) ;
241           if( sdigit->GetId() == curSDigit->GetId()) { // Are we in the same ECAL/HCAL/preshower tower ?              
242             *sdigit = *sdigit + *curSDigit;
243             newsdigit = kFALSE;
244           }
245         }
246         if (newsdigit) { 
247           new((*sdigits)[nSdigits])  AliEMCALDigit(*curSDigit);
248           nSdigits++ ;  
249         }
250         delete curSDigit ; 
251       }  // loop over all hits (hit = deposited energy/layer/entering particle)
252     } // loop over iprim
253     
254     sdigits->Sort() ;
255     
256     nSdigits = sdigits->GetEntriesFast() ;
257     fSDigitsInRun += nSdigits ;  
258     sdigits->Expand(nSdigits) ;
259
260     Int_t NPrimarymax = -1 ; 
261     Int_t i ;
262     for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) { 
263       AliEMCALDigit * sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(i)) ;
264       sdigit->SetIndexInList(i) ;
265     }
266     
267     for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) {   
268       if (((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) > NPrimarymax)
269         NPrimarymax = ((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) ;
270     }
271     
272     // Now write SDigits    
273     
274     //First list of sdigits
275
276     Int_t bufferSize = 32000 ;    
277     TBranch * sdigitsBranch = treeS->Branch("EMCAL",&sdigits,bufferSize);
278  
279     sdigitsBranch->Fill() ;
280
281     gime->WriteSDigits("OVERWRITE");
282     
283     //NEXT - SDigitizer
284
285     gime->WriteSDigitizer("OVERWRITE");
286     
287     if(strstr(option,"deb"))
288       PrintSDigits(option) ;  
289   }
290
291   Unload();
292   
293   gime->EmcalLoader()->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
294   
295   if(strstr(option,"tim")){
296     gBenchmark->Stop("EMCALSDigitizer"); 
297     Info("Exec", "took %f seconds for SDigitizing %f seconds per event", 
298          gBenchmark->GetCpuTime("EMCALSDigitizer"), gBenchmark->GetCpuTime("EMCALSDigitizer") ) ; 
299   }
300 }
301
302
303 //__________________________________________________________________
304 void AliEMCALSDigitizer::Print(Option_t*)const
305
306   // Prints parameters of SDigitizer
307   Info("Print", "\n------------------- %s -------------", GetName() ) ; 
308   printf("   Writing SDigits to branch with title  %s\n", fEventFolderName.Data()) ;
309   printf("   with digitization parameters  A = %f\n", fA) ; 
310   printf("                                 B = %f\n", fB) ;
311   printf("   Threshold for PRE Primary assignment= %f\n", fPREPrimThreshold)  ; 
312   printf("   Threshold for EC Primary assignment= %f\n", fECPrimThreshold)  ; 
313   printf("   Threshold for HC Primary assignment= %f\n", fHCPrimThreshold)  ; 
314   printf("---------------------------------------------------\n") ;
315
316 }
317
318 //__________________________________________________________________
319 Bool_t AliEMCALSDigitizer::operator==( AliEMCALSDigitizer const &sd )const
320 {
321   // Equal operator.
322   // SDititizers are equal if their pedestal, slope and threshold are equal
323
324   if( (fA==sd.fA)&&(fB==sd.fB)&&
325       (fECPrimThreshold==sd.fECPrimThreshold) &&
326       (fHCPrimThreshold==sd.fHCPrimThreshold) &&
327       (fPREPrimThreshold==sd.fPREPrimThreshold))
328     return kTRUE ;
329   else
330     return kFALSE ;
331 }
332
333 //__________________________________________________________________
334 void AliEMCALSDigitizer::PrintSDigits(Option_t * option){
335   //Prints list of digits produced at the current pass of AliEMCALDigitizer
336   
337   AliEMCALGetter * gime = AliEMCALGetter::Instance() ; 
338   const TClonesArray * sdigits = gime->SDigits() ; 
339   
340   TString message("\n") ;  
341   message += "event " ; 
342   message += gAlice->GetEvNumber() ;
343   message += "\n      Number of entries in SDigits list " ;
344   message += sdigits->GetEntriesFast() ; 
345   if(strstr(option,"all")||strstr(option,"EMC")){
346     
347     //loop over digits
348     AliEMCALDigit * digit;
349     message += "\n   Id  Amplitude    Time          Index Nprim: Primaries list \n" ;    
350     Int_t index ;
351     char * tempo = new char[8192]; 
352     for (index = 0 ; index < sdigits->GetEntries() ; index++) {
353       digit = dynamic_cast<AliEMCALDigit *>( sdigits->At(index) ) ;
354       sprintf(tempo, "\n%6d  %8d    %6.5e %4d      %2d :",
355               digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;  
356       message += tempo ; 
357       
358       Int_t iprimary;
359       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
360         sprintf(tempo, "%d ",digit->GetPrimary(iprimary+1) ) ; 
361         message += tempo ; 
362       }          
363     }
364     delete tempo ;
365   }
366   Info("PrintSDigits", message.Data() ) ; 
367 }
368
369 //____________________________________________________________________________ 
370 void AliEMCALSDigitizer::Unload() const
371 {
372   AliEMCALGetter * gime = AliEMCALGetter::Instance() ; 
373   AliEMCALLoader * loader = gime->EmcalLoader() ; 
374   loader->UnloadHits() ; 
375   loader->UnloadSDigits() ; 
376 }