]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALSDigitizer.cxx
Return to the original methods from AliLoader (Marco)
[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 <TH1.h>
54 #include <TBrowser.h>
55 #include <Riostream.h>
56 #include <TMath.h>
57
58 // --- Standard library ---
59 #include "stdlib.h"
60
61 // --- AliRoot header files ---
62 #include "AliLog.h"
63 #include "AliRunLoader.h"
64 #include "AliStack.h"
65 #include "AliEMCALDigit.h"
66 #include "AliEMCALLoader.h"
67 #include "AliEMCALHit.h"
68 #include "AliEMCALSDigitizer.h"
69 #include "AliEMCALGeometry.h"
70 #include "AliEMCALHistoUtilities.h"
71
72 ClassImp(AliEMCALSDigitizer)
73            
74 //____________________________________________________________________________ 
75   AliEMCALSDigitizer::AliEMCALSDigitizer():TTask("","") 
76 {
77   // ctor
78   fFirstEvent = fLastEvent  = fControlHists = 0;  
79   fDefaultInit = kTRUE ; 
80   fHists = 0;
81 }
82
83 //____________________________________________________________________________ 
84 AliEMCALSDigitizer::AliEMCALSDigitizer(const char * alirunFileName, 
85                                        const char * eventFolderName):
86   TTask("EMCAL"+AliConfig::Instance()->GetSDigitizerTaskName(), alirunFileName),
87   fEventFolderName(eventFolderName)
88 {
89   // ctor
90   fFirstEvent = fLastEvent  = fControlHists = 0 ; // runs one event by defaut  
91   Init();
92   InitParameters() ; 
93   fDefaultInit = kFALSE ; 
94   fHists = 0;
95   fControlHists = 1;
96   if(fControlHists) BookControlHists(1);
97 }
98
99
100 //____________________________________________________________________________ 
101 AliEMCALSDigitizer::AliEMCALSDigitizer(const AliEMCALSDigitizer & sd) : TTask(sd) {
102   //cpy ctor 
103
104   fFirstEvent    = sd.fFirstEvent ; 
105   fLastEvent     = sd.fLastEvent ;
106   fA             = sd.fA ;
107   fB             = sd.fB ;
108   fECPrimThreshold  = sd.fECPrimThreshold ;
109   fSDigitsInRun  = sd.fSDigitsInRun ;
110   SetName(sd.GetName()) ; 
111   SetTitle(sd.GetTitle()) ; 
112   fEventFolderName = sd.fEventFolderName;
113 }
114
115
116 //____________________________________________________________________________ 
117 AliEMCALSDigitizer::~AliEMCALSDigitizer() {
118   //dtor
119   AliLoader *emcalLoader=0;
120   if ((emcalLoader = AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL")))
121       emcalLoader->CleanSDigitizer();
122 }
123
124 //____________________________________________________________________________ 
125 void AliEMCALSDigitizer::Init(){
126   // Initialization: open root-file, allocate arrays for hits and sdigits,
127   // attach task SDigitizer to the list of EMCAL tasks
128   // 
129   // Initialization can not be done in the default constructor
130   //============================================================= YS
131   //  The initialisation is now done by the getter
132
133   fInit = kTRUE ; 
134    
135   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
136
137   if ( emcalLoader == 0 ) {
138     Fatal("Init", "Could not obtain the AliEMCALLoader");
139     return ;
140   } 
141   
142   emcalLoader->PostSDigitizer(this);
143   emcalLoader->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
144   
145 }
146
147 //____________________________________________________________________________ 
148 void AliEMCALSDigitizer::InitParameters()
149 {
150   //initialize parameters for sdigitization
151
152   const AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
153   if (geom->GetSampling() == 0.) {
154     Fatal("InitParameters", "Sampling factor not set !") ; 
155   }
156   // Initializes parameters
157   fA         = 0;
158   fB         = 1.e+7;
159   fSampling  = geom->GetSampling();
160
161  // threshold for deposit energy of hit
162   fECPrimThreshold  = 0.; // 24-nov-04 - was 1.e-6;
163   Print("");
164 }
165
166 //____________________________________________________________________________
167 void AliEMCALSDigitizer::Exec(Option_t *option) 
168
169   // Collects all hit of the same tower into digits
170   TString o(option); o.ToUpper();
171   if (strstr(option, "print") ) {
172     Print() ; 
173     return ; 
174   }
175   
176   if(strstr(option,"tim"))
177     gBenchmark->Start("EMCALSDigitizer");
178
179   AliRunLoader *rl = AliRunLoader::GetRunLoader();
180   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
181
182   //switch off reloading of this task while getting event
183   if (!fInit) { // to prevent overwrite existing file
184     AliError( Form("Give a version name different from %s", fEventFolderName.Data()) ) ;
185     return ;
186     }
187
188   if (fLastEvent == -1) 
189     fLastEvent = rl->GetNumberOfEvents() - 1 ;
190   else {
191     fLastEvent = TMath::Min(fLastEvent, rl->GetNumberOfEvents()-1);
192   }
193   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
194
195   Int_t ievent;
196   Float_t energy=0.; // de * fSampling - 23-nov-04
197   rl->LoadKinematics();
198   rl->LoadHits("EMCAL");
199
200   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
201     rl->GetEvent(ievent);
202     TTree * treeS = emcalLoader->TreeS();
203     if ( !treeS ) { 
204       emcalLoader->MakeSDigitsContainer();
205       treeS = emcalLoader->TreeS();
206     }
207     TClonesArray * hits = emcalLoader->Hits() ; 
208     TClonesArray * sdigits = emcalLoader->SDigits() ;
209     sdigits->Clear();
210     Int_t nSdigits = 0 ;
211     //Now make SDigits from hits, for EMCAL it is the same, so just copy    
212     Int_t nPrim =  static_cast<Int_t>((emcalLoader->TreeH())->GetEntries()) ; 
213     // This is not true: there is only one list of hits (MvL jan 2006)
214     // Attention nPrim is the number of primaries tracked by Geant 
215     // and this number could be different to the number of Primaries in TreeK; 
216     Int_t iprim ;
217     for ( iprim = 0 ; iprim < nPrim ; iprim++ ) { 
218       //=========== Get the EMCAL branch from Hits Tree for the Primary iprim
219       rl->Stack()->Particle(iprim) ;
220       Int_t i;
221       AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance(); 
222       for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
223         AliEMCALHit * hit = dynamic_cast<AliEMCALHit*>(hits->At(i)) ;
224         AliEMCALDigit * curSDigit = 0 ;
225         AliEMCALDigit * sdigit = 0 ;
226         Bool_t newsdigit = kTRUE; 
227
228         // hit->GetId() - Absolute Id number EMCAL segment
229         if(geom->CheckAbsCellId(hit->GetId())) { // was IsInECA(hit->GetId())
230           energy = hit->GetEnergy() * fSampling; // 23-nov-04
231           if(energy >  fECPrimThreshold )
232         // Assign primary number only if deposited energy is significant
233             curSDigit =  new AliEMCALDigit( hit->GetPrimary(),
234                                             hit->GetIparent(), hit->GetId(), 
235                                             Digitize(energy), hit->GetTime() ) ;
236           else
237             curSDigit =  new AliEMCALDigit( -1               , 
238                                             -1               ,
239                                             hit->GetId(), 
240                                             Digitize(energy), hit->GetTime() ) ;
241         } else {
242           Warning("Exec"," abs id %i is bad \n", hit->GetId());
243           newsdigit = kFALSE;
244           curSDigit = 0;
245         }
246
247         if(curSDigit != 0){
248           for(Int_t check= 0; check < nSdigits ; check++) {
249             sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(check)) ;
250             
251             if( sdigit->GetId() == curSDigit->GetId()) { // Are we in the same ECAL tower ?              
252               *sdigit = *sdigit + *curSDigit;
253               newsdigit = kFALSE;
254             }
255           }
256         }
257         if (newsdigit) {
258           new((*sdigits)[nSdigits])  AliEMCALDigit(*curSDigit);
259           nSdigits++ ;  
260         }
261         delete curSDigit ; 
262       }  // loop over all hits (hit = deposited energy/entering particle)
263     } // loop over iprim
264     sdigits->Sort() ;
265     
266     nSdigits = sdigits->GetEntriesFast() ;
267     fSDigitsInRun += nSdigits ;  
268
269     Int_t nPrimarymax = -1 ; 
270     Int_t i ;
271     Double_t e=0.,esum=0.;
272     AliEMCALHistoUtilities::FillH1(fHists, 0, double(sdigits->GetEntriesFast()));
273     for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) { 
274       AliEMCALDigit * sdigit = dynamic_cast<AliEMCALDigit *>(sdigits->At(i)) ;
275       sdigit->SetIndexInList(i) ;
276
277       AliEMCALHistoUtilities::FillH1(fHists, 2, double(sdigit->GetAmp()));
278       e = double(Calibrate(sdigit->GetAmp()));
279       esum += e;
280       AliEMCALHistoUtilities::FillH1(fHists, 3, e);
281       AliEMCALHistoUtilities::FillH1(fHists, 4, double(sdigit->GetId()));
282     }
283     if(esum>0.) AliEMCALHistoUtilities::FillH1(fHists, 1, esum);
284     
285     for (i = 0 ; i < sdigits->GetEntriesFast() ; i++) { // for what 
286       if (((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) > nPrimarymax)
287         nPrimarymax = ((dynamic_cast<AliEMCALDigit *>(sdigits->At(i)))->GetNprimary()) ;
288     }
289
290     // Now write SDigits    
291
292     Int_t bufferSize = 32000 ;    
293     TBranch * sdigitsBranch = treeS->GetBranch("EMCAL");
294     if (sdigitsBranch)
295       sdigitsBranch->SetAddress(&sdigits);
296     else
297       treeS->Branch("EMCAL",&sdigits,bufferSize);
298  
299     treeS->Fill();
300     
301     emcalLoader->WriteSDigits("OVERWRITE");
302     
303     //NEXT - SDigitizer
304     emcalLoader->WriteSDigitizer("OVERWRITE");  // why in event cycle ?
305     
306     if(strstr(option,"deb"))
307       PrintSDigits(option) ;  
308   }
309
310   Unload();
311   
312   emcalLoader->GetSDigitsDataLoader()->GetBaseTaskLoader()->SetDoNotReload(kTRUE);
313   
314   if(strstr(option,"tim")){
315     gBenchmark->Stop("EMCALSDigitizer"); 
316     printf("\n Exec: took %f seconds for SDigitizing %f seconds per event\n", 
317          gBenchmark->GetCpuTime("EMCALSDigitizer"), gBenchmark->GetCpuTime("EMCALSDigitizer")/nEvents ) ; 
318   }
319 }
320
321 //__________________________________________________________________
322
323 Int_t AliEMCALSDigitizer::Digitize(Float_t energy)const {
324   // Digitize the energy
325     Double_t aSignal = fA + energy*fB;
326     if (TMath::Abs(aSignal)>2147483647.0) { 
327       //PH 2147483647 is the max. integer
328       //PH This apparently is a problem which needs investigation
329       AliWarning(Form("Too big or too small energy %f",aSignal));
330       aSignal = TMath::Sign((Double_t)2147483647,aSignal);
331     }
332     return (Int_t ) aSignal;
333   }
334  
335
336 //__________________________________________________________________
337
338 void AliEMCALSDigitizer::Print1(Option_t * option)
339 {
340   Print(); 
341   PrintSDigits(option);
342 }
343
344 //__________________________________________________________________
345 void AliEMCALSDigitizer::Print(Option_t *option) const
346
347   // Prints parameters of SDigitizer
348   printf("Print: \n------------------- %s ------------- option %s\n", GetName() , option) ; 
349   printf("   fInit                                 %i\n", int(fInit));
350   printf("   fFirstEvent                           %i\n", fFirstEvent);
351   printf("   fLastEvent                            %i\n", fLastEvent);
352   printf("   Writing SDigits to branch with title  %s\n", fEventFolderName.Data()) ;
353   printf("   with digitization parameters       A = %f\n", fA) ; 
354   printf("                                      B = %f\n", fB) ;
355   printf("   Threshold for EC Primary assignment  = %f\n", fECPrimThreshold)  ;
356   printf("   Sampling                             = %f\n", fSampling);
357   printf("---------------------------------------------------\n") ;
358 }
359
360 //__________________________________________________________________
361 Bool_t AliEMCALSDigitizer::operator==( AliEMCALSDigitizer const &sd )const
362 {
363   // Equal operator.
364   // SDititizers are equal if their pedestal, slope and threshold are equal
365   if( (fA==sd.fA)&&(fB==sd.fB)&&
366       (fECPrimThreshold==sd.fECPrimThreshold))
367     return kTRUE ;
368   else
369     return kFALSE ;
370 }
371
372 //__________________________________________________________________
373 void AliEMCALSDigitizer::PrintSDigits(Option_t * option)
374 {
375   //Prints list of digits produced at the current pass of AliEMCALDigitizer
376   
377   
378   // AliEMCALGetter * gime = AliEMCALGetter::Instance() ; 
379   AliEMCALLoader *rl = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
380   const TClonesArray * sdigits = rl->SDigits() ; 
381   
382   printf("\n") ;  
383   printf("event %i", rl->GetRunLoader()->GetEventNumber());
384   printf(" Number of entries in SDigits list %i", sdigits->GetEntriesFast()); 
385   if(strstr(option,"all")||strstr(option,"EMC")){
386     
387     //loop over digits
388     AliEMCALDigit * digit;
389     printf("\n   Id  Amplitude    Time          Index Nprim: Primaries list \n") ;    
390     Int_t index, isum=0;
391     char * tempo = new char[8192]; 
392     for (index = 0 ; index < sdigits->GetEntries() ; index++) {
393       digit = dynamic_cast<AliEMCALDigit *>( sdigits->At(index) ) ;
394       sprintf(tempo, "\n%6d  %8d    %6.5e %4d      %2d :",
395               digit->GetId(), digit->GetAmp(), digit->GetTime(), digit->GetIndexInList(), digit->GetNprimary()) ;  
396       printf(tempo);
397       isum += digit->GetAmp();
398       
399       Int_t iprimary;
400       for (iprimary=0; iprimary<digit->GetNprimary(); iprimary++) {
401         sprintf(tempo, "%d ",digit->GetPrimary(iprimary+1) ) ; 
402         printf(tempo); 
403       }          
404     }
405     delete tempo ;
406     printf("\n** Sum %i : %10.3f GeV/c **\n ", isum, double(isum)*1.e-6);
407   } else printf("\n");
408 }
409
410 //____________________________________________________________________________ 
411 void AliEMCALSDigitizer::Unload() const
412 {
413   // Unload Hits and SDigits from the folder
414   AliEMCALLoader *rl = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
415   rl->UnloadHits() ; 
416   rl->UnloadSDigits() ; 
417 }
418
419 void AliEMCALSDigitizer::Browse(TBrowser* b)
420 {
421   if(fHists) b->Add(fHists);
422   TTask::Browse(b);
423 }
424
425 TList *AliEMCALSDigitizer::BookControlHists(int var)
426
427   //book histograms for monitoring sdigitization
428   // 22-nov-04
429   gROOT->cd();
430   const AliEMCALGeometry *geom = AliEMCALGeometry::GetInstance() ;
431   if(var>=1){
432     printf(" AliEMCALSDigitizer::BookControlHists() in action \n");
433     new TH1F("HSDigiN",  "#EMCAL  sdigits ", 1001, -0.5, 1000.5);
434     new TH1F("HSDigiSumEnergy","Sum.EMCAL energy", 1000, 0.0, 100.);
435     new TH1F("HSDigiAmp",  "EMCAL sdigits amplitude", 1000, 0., 2.e+9);
436     new TH1F("HSDigiEnergy","EMCAL cell energy", 1000, 0.0, 100.);
437     new TH1F("HSDigiAbsId","EMCAL absID for sdigits",
438     geom->GetNCells(), 0.5, Double_t(geom->GetNCells())+0.5);
439   }
440
441   fHists = AliEMCALHistoUtilities::MoveHistsToList("EmcalSDigiControlHists", kFALSE);
442   fHists = 0;
443
444   return fHists;
445 }
446
447 void AliEMCALSDigitizer::SaveHists(const char* name, Bool_t kSingleKey, const char* opt)
448 {
449   AliEMCALHistoUtilities::SaveListOfHists(fHists, name, kSingleKey, opt); 
450 }