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